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-2011, 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 <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;
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: PetscTypeCompare((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: if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
211: IPView(qep->ip,viewer);
212: return(0);
213: }
217: /*@
218: QEPPrintSolution - Prints the computed eigenvalues.
220: Collective on QEP
222: Input Parameters:
223: + qep - the eigensolver context
224: - viewer - optional visualization context
226: Options Database:
227: . -qep_terse - print only minimal information
229: Note:
230: By default, this function prints a table with eigenvalues and associated
231: relative errors. With -qep_terse only the eigenvalues are printed.
233: Level: intermediate
235: .seealso: PetscViewerASCIIOpen()
236: @*/
237: PetscErrorCode QEPPrintSolution(QEP qep,PetscViewer viewer)
238: {
239: PetscBool terse,errok,isascii;
240: PetscReal error,re,im;
241: PetscScalar kr,ki;
242: PetscInt i,j;
247: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)qep)->comm);
250: if (!qep->eigr || !qep->eigi || !qep->V) {
251: SETERRQ(((PetscObject)qep)->comm,PETSC_ERR_ARG_WRONGSTATE,"QEPSolve must be called first");
252: }
253: PetscTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&isascii);
254: if (!isascii) return(0);
256: PetscOptionsHasName(PETSC_NULL,"-qep_terse",&terse);
257: if (terse) {
258: if (qep->nconv<qep->nev) {
259: PetscViewerASCIIPrintf(viewer," Problem: less than %D eigenvalues converged\n\n",qep->nev);
260: } else {
261: errok = PETSC_TRUE;
262: for (i=0;i<qep->nev;i++) {
263: QEPComputeRelativeError(qep,i,&error);
264: errok = (errok && error<qep->tol)? PETSC_TRUE: PETSC_FALSE;
265: }
266: if (errok) {
267: PetscViewerASCIIPrintf(viewer," All requested eigenvalues computed up to the required tolerance:");
268: for (i=0;i<=(qep->nev-1)/8;i++) {
269: PetscViewerASCIIPrintf(viewer,"\n ");
270: for (j=0;j<PetscMin(8,qep->nev-8*i);j++) {
271: QEPGetEigenpair(qep,8*i+j,&kr,&ki,PETSC_NULL,PETSC_NULL);
272: #if defined(PETSC_USE_COMPLEX)
273: re = PetscRealPart(kr);
274: im = PetscImaginaryPart(kr);
275: #else
276: re = kr;
277: im = ki;
278: #endif
279: if (PetscAbs(re)/PetscAbs(im)<PETSC_SMALL) re = 0.0;
280: if (PetscAbs(im)/PetscAbs(re)<PETSC_SMALL) im = 0.0;
281: if (im!=0.0) {
282: PetscViewerASCIIPrintf(viewer,"%.5F%+.5Fi",re,im);
283: } else {
284: PetscViewerASCIIPrintf(viewer,"%.5F",re);
285: }
286: if (8*i+j+1<qep->nev) { PetscViewerASCIIPrintf(viewer,", "); }
287: }
288: }
289: PetscViewerASCIIPrintf(viewer,"\n\n");
290: } else {
291: PetscViewerASCIIPrintf(viewer," Problem: some of the first %D relative errors are higher than the tolerance\n\n",qep->nev);
292: }
293: }
294: } else {
295: PetscViewerASCIIPrintf(viewer," Number of converged approximate eigenpairs: %D\n\n",qep->nconv);
296: if (qep->nconv>0) {
297: PetscViewerASCIIPrintf(viewer,
298: " k ||(k^2M+Ck+K)x||/||kx||\n"
299: " ----------------- -------------------------\n");
300: for (i=0;i<qep->nconv;i++) {
301: QEPGetEigenpair(qep,i,&kr,&ki,PETSC_NULL,PETSC_NULL);
302: QEPComputeRelativeError(qep,i,&error);
303: #if defined(PETSC_USE_COMPLEX)
304: re = PetscRealPart(kr);
305: im = PetscImaginaryPart(kr);
306: #else
307: re = kr;
308: im = ki;
309: #endif
310: if (im!=0.0) {
311: PetscViewerASCIIPrintf(viewer," % 9F%+9F i %12G\n",re,im,error);
312: } else {
313: PetscViewerASCIIPrintf(viewer," % 12F %12G\n",re,error);
314: }
315: }
316: PetscViewerASCIIPrintf(viewer,"\n");
317: }
318: }
319: return(0);
320: }
324: /*@C
325: QEPCreate - Creates the default QEP context.
327: Collective on MPI_Comm
329: Input Parameter:
330: . comm - MPI communicator
332: Output Parameter:
333: . qep - location to put the QEP context
335: Note:
336: The default QEP type is QEPLINEAR
338: Level: beginner
340: .seealso: QEPSetUp(), QEPSolve(), QEPDestroy(), QEP
341: @*/
342: PetscErrorCode QEPCreate(MPI_Comm comm,QEP *outqep)
343: {
345: QEP qep;
349: *outqep = 0;
350: PetscHeaderCreate(qep,_p_QEP,struct _QEPOps,QEP_CLASSID,-1,"QEP","Quadratic Eigenvalue Problem","QEP",comm,QEPDestroy,QEPView);
352: qep->M = 0;
353: qep->C = 0;
354: qep->K = 0;
355: qep->max_it = 0;
356: qep->nev = 1;
357: qep->ncv = 0;
358: qep->mpd = 0;
359: qep->nini = 0;
360: qep->ninil = 0;
361: qep->allocated_ncv = 0;
362: qep->tol = PETSC_DEFAULT;
363: qep->sfactor = 0.0;
364: qep->conv_func = QEPDefaultConverged;
365: qep->conv_ctx = PETSC_NULL;
366: qep->which = (QEPWhich)0;
367: qep->which_func = PETSC_NULL;
368: qep->which_ctx = PETSC_NULL;
369: qep->leftvecs = PETSC_FALSE;
370: qep->problem_type = (QEPProblemType)0;
371: qep->V = PETSC_NULL;
372: qep->W = PETSC_NULL;
373: qep->IS = PETSC_NULL;
374: qep->ISL = PETSC_NULL;
375: qep->T = PETSC_NULL;
376: qep->eigr = PETSC_NULL;
377: qep->eigi = PETSC_NULL;
378: qep->errest = PETSC_NULL;
379: qep->data = PETSC_NULL;
380: qep->t = PETSC_NULL;
381: qep->nconv = 0;
382: qep->its = 0;
383: qep->perm = PETSC_NULL;
384: qep->matvecs = 0;
385: qep->linits = 0;
386: qep->nwork = 0;
387: qep->work = PETSC_NULL;
388: qep->setupcalled = 0;
389: qep->reason = QEP_CONVERGED_ITERATING;
390: qep->numbermonitors = 0;
391: qep->trackall = PETSC_FALSE;
392: qep->rand = 0;
394: PetscRandomCreate(comm,&qep->rand);
395: PetscRandomSetSeed(qep->rand,0x12345678);
396: PetscLogObjectParent(qep,qep->rand);
397: *outqep = qep;
398: return(0);
399: }
400:
403: /*@C
404: QEPSetType - Selects the particular solver to be used in the QEP object.
406: Logically Collective on QEP
408: Input Parameters:
409: + qep - the quadratic eigensolver context
410: - type - a known method
412: Options Database Key:
413: . -qep_type <method> - Sets the method; use -help for a list
414: of available methods
415:
416: Notes:
417: See "slepc/include/slepcqep.h" for available methods. The default
418: is QEPLINEAR.
420: Normally, it is best to use the QEPSetFromOptions() command and
421: then set the QEP type from the options database rather than by using
422: this routine. Using the options database provides the user with
423: maximum flexibility in evaluating the different available methods.
424: The QEPSetType() routine is provided for those situations where it
425: is necessary to set the iterative solver independently of the command
426: line or options database.
428: Level: intermediate
430: .seealso: QEPType
431: @*/
432: PetscErrorCode QEPSetType(QEP qep,const QEPType type)
433: {
434: PetscErrorCode ierr,(*r)(QEP);
435: PetscBool match;
441: PetscTypeCompare((PetscObject)qep,type,&match);
442: if (match) return(0);
444: PetscFListFind(QEPList,((PetscObject)qep)->comm,type,PETSC_TRUE,(void (**)(void))&r);
445: if (!r) SETERRQ1(((PetscObject)qep)->comm,PETSC_ERR_ARG_UNKNOWN_TYPE,"Unknown QEP type given: %s",type);
447: if (qep->ops->destroy) { (*qep->ops->destroy)(qep); }
448: PetscMemzero(qep->ops,sizeof(struct _QEPOps));
450: qep->setupcalled = 0;
451: PetscObjectChangeTypeName((PetscObject)qep,type);
452: (*r)(qep);
453: return(0);
454: }
458: /*@C
459: QEPGetType - Gets the QEP type as a string from the QEP object.
461: Not Collective
463: Input Parameter:
464: . qep - the eigensolver context
466: Output Parameter:
467: . name - name of QEP method
469: Level: intermediate
471: .seealso: QEPSetType()
472: @*/
473: PetscErrorCode QEPGetType(QEP qep,const QEPType *type)
474: {
478: *type = ((PetscObject)qep)->type_name;
479: return(0);
480: }
484: /*@C
485: QEPRegister - See QEPRegisterDynamic()
487: Level: advanced
488: @*/
489: PetscErrorCode QEPRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(QEP))
490: {
492: char fullname[PETSC_MAX_PATH_LEN];
495: PetscFListConcat(path,name,fullname);
496: PetscFListAdd(&QEPList,sname,fullname,(void (*)(void))function);
497: return(0);
498: }
502: /*@
503: QEPRegisterDestroy - Frees the list of QEP methods that were
504: registered by QEPRegisterDynamic().
506: Not Collective
508: Level: advanced
510: .seealso: QEPRegisterDynamic(), QEPRegisterAll()
511: @*/
512: PetscErrorCode QEPRegisterDestroy(void)
513: {
517: PetscFListDestroy(&QEPList);
518: QEPRegisterAllCalled = PETSC_FALSE;
519: return(0);
520: }
524: /*@C
525: QEPReset - Resets the QEP context to the setupcalled=0 state and removes any
526: allocated objects.
528: Collective on QEP
530: Input Parameter:
531: . qep - eigensolver context obtained from QEPCreate()
533: Level: advanced
535: .seealso: QEPDestroy()
536: @*/
537: PetscErrorCode QEPReset(QEP qep)
538: {
543: if (qep->ops->reset) { (qep->ops->reset)(qep); }
544: if (qep->ip) { IPReset(qep->ip); }
545: MatDestroy(&qep->M);
546: MatDestroy(&qep->C);
547: MatDestroy(&qep->K);
548: VecDestroy(&qep->t);
549: QEPFreeSolution(qep);
550: qep->matvecs = 0;
551: qep->linits = 0;
552: qep->setupcalled = 0;
553: return(0);
554: }
558: /*@C
559: QEPDestroy - Destroys the QEP context.
561: Collective on QEP
563: Input Parameter:
564: . qep - eigensolver context obtained from QEPCreate()
566: Level: beginner
568: .seealso: QEPCreate(), QEPSetUp(), QEPSolve()
569: @*/
570: PetscErrorCode QEPDestroy(QEP *qep)
571: {
575: if (!*qep) return(0);
577: if (--((PetscObject)(*qep))->refct > 0) { *qep = 0; return(0); }
578: QEPReset(*qep);
579: PetscObjectDepublish(*qep);
580: if ((*qep)->ops->destroy) { (*(*qep)->ops->destroy)(*qep); }
581: IPDestroy(&(*qep)->ip);
582: PetscRandomDestroy(&(*qep)->rand);
583: QEPMonitorCancel(*qep);
584: PetscHeaderDestroy(qep);
585: return(0);
586: }
590: /*@
591: QEPSetIP - Associates an inner product object to the quadratic eigensolver.
593: Collective on QEP
595: Input Parameters:
596: + qep - eigensolver context obtained from QEPCreate()
597: - ip - the inner product object
599: Note:
600: Use QEPGetIP() to retrieve the inner product context (for example,
601: to free it at the end of the computations).
603: Level: advanced
605: .seealso: QEPGetIP()
606: @*/
607: PetscErrorCode QEPSetIP(QEP qep,IP ip)
608: {
615: PetscObjectReference((PetscObject)ip);
616: IPDestroy(&qep->ip);
617: qep->ip = ip;
618: PetscLogObjectParent(qep,qep->ip);
619: return(0);
620: }
624: /*@C
625: QEPGetIP - Obtain the inner product object associated
626: to the quadratic eigensolver object.
628: Not Collective
630: Input Parameters:
631: . qep - eigensolver context obtained from QEPCreate()
633: Output Parameter:
634: . ip - inner product context
636: Level: advanced
638: .seealso: QEPSetIP()
639: @*/
640: PetscErrorCode QEPGetIP(QEP qep,IP *ip)
641: {
647: if (!qep->ip) {
648: IPCreate(((PetscObject)qep)->comm,&qep->ip);
649: PetscLogObjectParent(qep,qep->ip);
650: }
651: *ip = qep->ip;
652: return(0);
653: }