Actual source code: qepbasic.c
1: /*
2: The basic QEP routines, Create, View, etc. are here.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2012, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
9:
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc-private/qepimpl.h> /*I "slepcqep.h" I*/
26: PetscFList QEPList = 0;
27: PetscBool QEPRegisterAllCalled = PETSC_FALSE;
28: PetscClassId QEP_CLASSID = 0;
29: PetscLogEvent QEP_SetUp = 0,QEP_Solve = 0,QEP_Dense = 0;
30: static PetscBool QEPPackageInitialized = PETSC_FALSE;
34: /*@C
35: QEPFinalizePackage - This function destroys everything in the Slepc interface
36: to the QEP package. It is called from SlepcFinalize().
38: Level: developer
40: .seealso: SlepcFinalize()
41: @*/
42: PetscErrorCode QEPFinalizePackage(void)
43: {
45: QEPPackageInitialized = PETSC_FALSE;
46: QEPList = 0;
47: QEPRegisterAllCalled = PETSC_FALSE;
48: return(0);
49: }
53: /*@C
54: QEPInitializePackage - This function initializes everything in the QEP package. It is called
55: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to QEPCreate()
56: when using static libraries.
58: Input Parameter:
59: . path - The dynamic library path, or PETSC_NULL
61: Level: developer
63: .seealso: SlepcInitialize()
64: @*/
65: PetscErrorCode QEPInitializePackage(const char *path)
66: {
67: char logList[256];
68: char *className;
69: PetscBool opt;
73: if (QEPPackageInitialized) return(0);
74: QEPPackageInitialized = PETSC_TRUE;
75: /* Register Classes */
76: PetscClassIdRegister("Quadratic Eigenproblem Solver",&QEP_CLASSID);
77: /* Register Constructors */
78: QEPRegisterAll(path);
79: /* Register Events */
80: PetscLogEventRegister("QEPSetUp",QEP_CLASSID,&QEP_SetUp);
81: PetscLogEventRegister("QEPSolve",QEP_CLASSID,&QEP_Solve);
82: PetscLogEventRegister("QEPDense",QEP_CLASSID,&QEP_Dense);
83: /* Process info exclusions */
84: PetscOptionsGetString(PETSC_NULL,"-info_exclude",logList,256,&opt);
85: if (opt) {
86: PetscStrstr(logList,"qep",&className);
87: if (className) {
88: PetscInfoDeactivateClass(QEP_CLASSID);
89: }
90: }
91: /* Process summary exclusions */
92: PetscOptionsGetString(PETSC_NULL,"-log_summary_exclude",logList,256,&opt);
93: if (opt) {
94: PetscStrstr(logList,"qep",&className);
95: if (className) {
96: PetscLogEventDeactivateClass(QEP_CLASSID);
97: }
98: }
99: PetscRegisterFinalize(QEPFinalizePackage);
100: return(0);
101: }
105: /*@C
106: QEPView - Prints the QEP data structure.
108: Collective on QEP
110: Input Parameters:
111: + qep - the quadratic eigenproblem solver context
112: - viewer - optional visualization context
114: Options Database Key:
115: . -qep_view - Calls QEPView() at end of QEPSolve()
117: Note:
118: The available visualization contexts include
119: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
120: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
121: output where only the first processor opens
122: the file. All other processors send their
123: data to the first processor to print.
125: The user can open an alternative visualization context with
126: PetscViewerASCIIOpen() - output to a specified file.
128: Level: beginner
130: .seealso: PetscViewerASCIIOpen()
131: @*/
132: PetscErrorCode QEPView(QEP qep,PetscViewer viewer)
133: {
135: const char *type;
136: PetscBool isascii,islinear;
140: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)qep)->comm);
144: #if defined(PETSC_USE_COMPLEX)
145: #define HERM "hermitian"
146: #else
147: #define HERM "symmetric"
148: #endif
149: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
150: if (isascii) {
151: PetscObjectPrintClassNamePrefixType((PetscObject)qep,viewer,"QEP Object");
152: if (qep->ops->view) {
153: PetscViewerASCIIPushTab(viewer);
154: (*qep->ops->view)(qep,viewer);
155: PetscViewerASCIIPopTab(viewer);
156: }
157: if (qep->problem_type) {
158: switch (qep->problem_type) {
159: case QEP_GENERAL: type = "general quadratic eigenvalue problem"; break;
160: case QEP_HERMITIAN: type = HERM " quadratic eigenvalue problem"; break;
161: case QEP_GYROSCOPIC: type = "gyroscopic quadratic eigenvalue problem"; break;
162: default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->problem_type");
163: }
164: } else type = "not yet set";
165: PetscViewerASCIIPrintf(viewer," problem type: %s\n",type);
166: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: ");
167: if (!qep->which) {
168: PetscViewerASCIIPrintf(viewer,"not yet set\n");
169: } else switch (qep->which) {
170: case QEP_LARGEST_MAGNITUDE:
171: PetscViewerASCIIPrintf(viewer,"largest eigenvalues in magnitude\n");
172: break;
173: case QEP_SMALLEST_MAGNITUDE:
174: PetscViewerASCIIPrintf(viewer,"smallest eigenvalues in magnitude\n");
175: break;
176: case QEP_LARGEST_REAL:
177: PetscViewerASCIIPrintf(viewer,"largest real parts\n");
178: break;
179: case QEP_SMALLEST_REAL:
180: PetscViewerASCIIPrintf(viewer,"smallest real parts\n");
181: break;
182: case QEP_LARGEST_IMAGINARY:
183: PetscViewerASCIIPrintf(viewer,"largest imaginary parts\n");
184: break;
185: case QEP_SMALLEST_IMAGINARY:
186: PetscViewerASCIIPrintf(viewer,"smallest imaginary parts\n");
187: break;
188: default: SETERRQ(((PetscObject)qep)->comm,1,"Wrong value of qep->which");
189: }
190: if (qep->leftvecs) {
191: PetscViewerASCIIPrintf(viewer," computing left eigenvectors also\n");
192: }
193: PetscViewerASCIIPrintf(viewer," number of eigenvalues (nev): %D\n",qep->nev);
194: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %D\n",qep->ncv);
195: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %D\n",qep->mpd);
196: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %D\n",qep->max_it);
197: PetscViewerASCIIPrintf(viewer," tolerance: %G\n",qep->tol);
198: PetscViewerASCIIPrintf(viewer," scaling factor: %G\n",qep->sfactor);
199: if (qep->nini!=0) {
200: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial space: %D\n",PetscAbs(qep->nini));
201: }
202: if (qep->ninil!=0) {
203: PetscViewerASCIIPrintf(viewer," dimension of user-provided initial left space: %D\n",PetscAbs(qep->ninil));
204: }
205: } else {
206: if (qep->ops->view) {
207: (*qep->ops->view)(qep,viewer);
208: }
209: }
210: PetscObjectTypeCompare((PetscObject)qep,QEPLINEAR,&islinear);
211: if (!islinear) {
212: if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
213: IPView(qep->ip,viewer);
214: if (!qep->ds) { QEPGetDS(qep,&qep->ds); }
215: PetscViewerPushFormat(viewer,PETSC_VIEWER_ASCII_INFO);
216: DSView(qep->ds,viewer);
217: PetscViewerPopFormat(viewer);
218: }
219: return(0);
220: }
224: /*@
225: QEPPrintSolution - Prints the computed eigenvalues.
227: Collective on QEP
229: Input Parameters:
230: + qep - the eigensolver context
231: - viewer - optional visualization context
233: Options Database:
234: . -qep_terse - print only minimal information
236: Note:
237: By default, this function prints a table with eigenvalues and associated
238: relative errors. With -qep_terse only the eigenvalues are printed.
240: Level: intermediate
242: .seealso: PetscViewerASCIIOpen()
243: @*/
244: PetscErrorCode QEPPrintSolution(QEP qep,PetscViewer viewer)
245: {
246: PetscBool terse,errok,isascii;
247: PetscReal error,re,im;
248: PetscScalar kr,ki;
249: PetscInt i,j;
254: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)qep)->comm);
257: if (!qep->eigr || !qep->eigi || !qep->V) SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSolve must be called first");
258: PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
259: if (!isascii) return(0);
261: PetscOptionsHasName(PETSC_NULL,"-qep_terse",&terse);
262: if (terse) {
263: if (qep->nconv<qep->nev) {
264: PetscViewerASCIIPrintf(viewer," Problem: less than %D eigenvalues converged\n\n",qep->nev);
265: } else {
266: errok = PETSC_TRUE;
267: for (i=0;i<qep->nev;i++) {
268: QEPComputeRelativeError(qep,i,&error);
269: errok = (errok && error<qep->tol)? PETSC_TRUE: PETSC_FALSE;
270: }
271: if (errok) {
272: PetscViewerASCIIPrintf(viewer," All requested eigenvalues computed up to the required tolerance:");
273: for (i=0;i<=(qep->nev-1)/8;i++) {
274: PetscViewerASCIIPrintf(viewer,"\n ");
275: for (j=0;j<PetscMin(8,qep->nev-8*i);j++) {
276: QEPGetEigenpair(qep,8*i+j,&kr,&ki,PETSC_NULL,PETSC_NULL);
277: #if defined(PETSC_USE_COMPLEX)
278: re = PetscRealPart(kr);
279: im = PetscImaginaryPart(kr);
280: #else
281: re = kr;
282: im = ki;
283: #endif
284: if (PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
285: if (PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
286: if (im!=0.0) {
287: PetscViewerASCIIPrintf(viewer,"%.5F%+.5Fi",re,im);
288: } else {
289: PetscViewerASCIIPrintf(viewer,"%.5F",re);
290: }
291: if (8*i+j+1<qep->nev) { PetscViewerASCIIPrintf(viewer,", "); }
292: }
293: }
294: PetscViewerASCIIPrintf(viewer,"\n\n");
295: } else {
296: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",qep->nev);
297: }
298: }
299: } else {
300: PetscViewerASCIIPrintf(viewer," Number of converged approximate eigenpairs: %D\n\n",qep->nconv);
301: if (qep->nconv>0) {
302: PetscViewerASCIIPrintf(viewer,
303: " k ||(k^2M+Ck+K)x||/||kx||\n"
304: " ----------------- -------------------------\n");
305: for (i=0;i<qep->nconv;i++) {
306: QEPGetEigenpair(qep,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
307: QEPComputeRelativeError(qep,i,&error);
308: #if defined(PETSC_USE_COMPLEX)
309: re = PetscRealPart(kr);
310: im = PetscImaginaryPart(kr);
311: #else
312: re = kr;
313: im = ki;
314: #endif
315: if (im!=0.0) {
316: PetscViewerASCIIPrintf(viewer," % 9F%+9F i %12G\n",re,im,error);
317: } else {
318: PetscViewerASCIIPrintf(viewer," % 12F %12G\n",re,error);
319: }
320: }
321: PetscViewerASCIIPrintf(viewer,"\n");
322: }
323: }
324: return(0);
325: }
329: /*@C
330: QEPCreate - Creates the default QEP context.
332: Collective on MPI_Comm
334: Input Parameter:
335: . comm - MPI communicator
337: Output Parameter:
338: . qep - location to put the QEP context
340: Note:
341: The default QEP type is QEPLINEAR
343: Level: beginner
345: .seealso: QEPSetUp(), QEPSolve(), QEPDestroy(), QEP
346: @*/
347: PetscErrorCode QEPCreate(MPI_Comm comm,QEP *outqep)
348: {
350: QEP qep;
354: *outqep = 0;
355: SlepcHeaderCreate(qep,_p_QEP,struct _QEPOps,QEP_CLASSID,-1,"QEP","Quadratic Eigenvalue Problem","QEP",comm,QEPDestroy,QEPView);
357: qep->M = 0;
358: qep->C = 0;
359: qep->K = 0;
360: qep->max_it = 0;
361: qep->nev = 1;
362: qep->ncv = 0;
363: qep->mpd = 0;
364: qep->nini = 0;
365: qep->ninil = 0;
366: qep->allocated_ncv = 0;
367: qep->ip = 0;
368: qep->ds = 0;
369: qep->tol = PETSC_DEFAULT;
370: qep->sfactor = 0.0;
371: qep->conv_func = QEPConvergedDefault;
372: qep->conv_ctx = PETSC_NULL;
373: qep->which = (QEPWhich)0;
374: qep->which_func = PETSC_NULL;
375: qep->which_ctx = PETSC_NULL;
376: qep->leftvecs = PETSC_FALSE;
377: qep->problem_type = (QEPProblemType)0;
378: qep->V = PETSC_NULL;
379: qep->W = PETSC_NULL;
380: qep->IS = PETSC_NULL;
381: qep->ISL = PETSC_NULL;
382: qep->eigr = PETSC_NULL;
383: qep->eigi = PETSC_NULL;
384: qep->errest = PETSC_NULL;
385: qep->data = PETSC_NULL;
386: qep->t = PETSC_NULL;
387: qep->nconv = 0;
388: qep->its = 0;
389: qep->perm = PETSC_NULL;
390: qep->matvecs = 0;
391: qep->linits = 0;
392: qep->nwork = 0;
393: qep->work = PETSC_NULL;
394: qep->setupcalled = 0;
395: qep->reason = QEP_CONVERGED_ITERATING;
396: qep->numbermonitors = 0;
397: qep->trackall = PETSC_FALSE;
398: qep->rand = 0;
400: PetscRandomCreate(comm,&qep->rand);
401: PetscRandomSetSeed(qep->rand,0x12345678);
402: PetscLogObjectParent(qep,qep->rand);
403: *outqep = qep;
404: return(0);
405: }
406:
409: /*@C
410: QEPSetType - Selects the particular solver to be used in the QEP object.
412: Logically Collective on QEP
414: Input Parameters:
415: + qep - the quadratic eigensolver context
416: - type - a known method
418: Options Database Key:
419: . -qep_type <method> - Sets the method; use -help for a list
420: of available methods
421:
422: Notes:
423: See "slepc/include/slepcqep.h" for available methods. The default
424: is QEPLINEAR.
426: Normally, it is best to use the QEPSetFromOptions() command and
427: then set the QEP type from the options database rather than by using
428: this routine. Using the options database provides the user with
429: maximum flexibility in evaluating the different available methods.
430: The QEPSetType() routine is provided for those situations where it
431: is necessary to set the iterative solver independently of the command
432: line or options database.
434: Level: intermediate
436: .seealso: QEPType
437: @*/
438: PetscErrorCode QEPSetType(QEP qep,const QEPType type)
439: {
440: PetscErrorCode ierr,(*r)(QEP);
441: PetscBool match;
447: PetscObjectTypeCompare((PetscObject)qep,type,&match);
448: if (match) return(0);
450: PetscFListFind(QEPList,((PetscObject)qep)->comm,type,PETSC_TRUE,(void (**)(void))&r);
451: if (!r) SETERRQ1(((PetscObject)qep)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown QEP type given: %s",type);
453: if (qep->ops->destroy) { (*qep->ops->destroy)(qep); }
454: PetscMemzero(qep->ops,sizeof(struct _QEPOps));
456: qep->setupcalled = 0;
457: PetscObjectChangeTypeName((PetscObject)qep,type);
458: (*r)(qep);
459: return(0);
460: }
464: /*@C
465: QEPGetType - Gets the QEP type as a string from the QEP object.
467: Not Collective
469: Input Parameter:
470: . qep - the eigensolver context
472: Output Parameter:
473: . name - name of QEP method
475: Level: intermediate
477: .seealso: QEPSetType()
478: @*/
479: PetscErrorCode QEPGetType(QEP qep,const QEPType *type)
480: {
484: *type = ((PetscObject)qep)->type_name;
485: return(0);
486: }
490: /*@C
491: QEPRegister - See QEPRegisterDynamic()
493: Level: advanced
494: @*/
495: PetscErrorCode QEPRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(QEP))
496: {
498: char fullname[PETSC_MAX_PATH_LEN];
501: PetscFListConcat(path,name,fullname);
502: PetscFListAdd(&QEPList,sname,fullname,(void (*)(void))function);
503: return(0);
504: }
508: /*@
509: QEPRegisterDestroy - Frees the list of QEP methods that were
510: registered by QEPRegisterDynamic().
512: Not Collective
514: Level: advanced
516: .seealso: QEPRegisterDynamic(), QEPRegisterAll()
517: @*/
518: PetscErrorCode QEPRegisterDestroy(void)
519: {
523: PetscFListDestroy(&QEPList);
524: QEPRegisterAllCalled = PETSC_FALSE;
525: return(0);
526: }
530: /*@
531: QEPReset - Resets the QEP context to the setupcalled=0 state and removes any
532: allocated objects.
534: Collective on QEP
536: Input Parameter:
537: . qep - eigensolver context obtained from QEPCreate()
539: Level: advanced
541: .seealso: QEPDestroy()
542: @*/
543: PetscErrorCode QEPReset(QEP qep)
544: {
549: if (qep->ops->reset) { (qep->ops->reset)(qep); }
550: if (qep->ip) { IPReset(qep->ip); }
551: if (qep->ds) { DSReset(qep->ds); }
552: MatDestroy(&qep->M);
553: MatDestroy(&qep->C);
554: MatDestroy(&qep->K);
555: VecDestroy(&qep->t);
556: QEPFreeSolution(qep);
557: qep->matvecs = 0;
558: qep->linits = 0;
559: qep->setupcalled = 0;
560: return(0);
561: }
565: /*@C
566: QEPDestroy - Destroys the QEP context.
568: Collective on QEP
570: Input Parameter:
571: . qep - eigensolver context obtained from QEPCreate()
573: Level: beginner
575: .seealso: QEPCreate(), QEPSetUp(), QEPSolve()
576: @*/
577: PetscErrorCode QEPDestroy(QEP *qep)
578: {
582: if (!*qep) return(0);
584: if (--((PetscObject)(*qep))->refct > 0) { *qep = 0; return(0); }
585: QEPReset(*qep);
586: PetscObjectDepublish(*qep);
587: if ((*qep)->ops->destroy) { (*(*qep)->ops->destroy)(*qep); }
588: IPDestroy(&(*qep)->ip);
589: DSDestroy(&(*qep)->ds);
590: PetscRandomDestroy(&(*qep)->rand);
591: QEPMonitorCancel(*qep);
592: PetscHeaderDestroy(qep);
593: return(0);
594: }
598: /*@
599: QEPSetIP - Associates an inner product object to the quadratic eigensolver.
601: Collective on QEP
603: Input Parameters:
604: + qep - eigensolver context obtained from QEPCreate()
605: - ip - the inner product object
607: Note:
608: Use QEPGetIP() to retrieve the inner product context (for example,
609: to free it at the end of the computations).
611: Level: advanced
613: .seealso: QEPGetIP()
614: @*/
615: PetscErrorCode QEPSetIP(QEP qep,IP ip)
616: {
623: PetscObjectReference((PetscObject)ip);
624: IPDestroy(&qep->ip);
625: qep->ip = ip;
626: PetscLogObjectParent(qep,qep->ip);
627: return(0);
628: }
632: /*@C
633: QEPGetIP - Obtain the inner product object associated
634: to the quadratic eigensolver object.
636: Not Collective
638: Input Parameters:
639: . qep - eigensolver context obtained from QEPCreate()
641: Output Parameter:
642: . ip - inner product context
644: Level: advanced
646: .seealso: QEPSetIP()
647: @*/
648: PetscErrorCode QEPGetIP(QEP qep,IP *ip)
649: {
655: if (!qep->ip) {
656: IPCreate(((PetscObject)qep)->comm,&qep->ip);
657: PetscLogObjectParent(qep,qep->ip);
658: }
659: *ip = qep->ip;
660: return(0);
661: }
665: /*@
666: QEPSetDS - Associates a direct solver object to the quadratic eigensolver.
668: Collective on QEP
670: Input Parameters:
671: + qep - eigensolver context obtained from QEPCreate()
672: - ds - the direct solver object
674: Note:
675: Use QEPGetDS() to retrieve the direct solver context (for example,
676: to free it at the end of the computations).
678: Level: advanced
680: .seealso: QEPGetDS()
681: @*/
682: PetscErrorCode QEPSetDS(QEP qep,DS ds)
683: {
690: PetscObjectReference((PetscObject)ds);
691: DSDestroy(&qep->ds);
692: qep->ds = ds;
693: PetscLogObjectParent(qep,qep->ds);
694: return(0);
695: }
699: /*@C
700: QEPGetDS - Obtain the direct solver object associated to the
701: quadratic eigensolver object.
703: Not Collective
705: Input Parameters:
706: . qep - eigensolver context obtained from QEPCreate()
708: Output Parameter:
709: . ds - direct solver context
711: Level: advanced
713: .seealso: QEPSetDS()
714: @*/
715: PetscErrorCode QEPGetDS(QEP qep,DS *ds)
716: {
722: if (!qep->ds) {
723: DSCreate(((PetscObject)qep)->comm,&qep->ds);
724: PetscLogObjectParent(qep,qep->ds);
725: }
726: *ds = qep->ds;
727: return(0);
728: }