Actual source code: basic.c
1: /*
2: The basic EPS routines, Create, View, etc. are here.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2009, Universidad 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/epsimpl.h
26: PetscFList EPSList = 0;
27: PetscCookie EPS_COOKIE = 0;
28: PetscLogEvent EPS_SetUp = 0, EPS_Solve = 0, EPS_Dense = 0;
32: /*@C
33: EPSInitializePackage - This function initializes everything in the EPS package. It is called
34: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to EPSCreate()
35: when using static libraries.
37: Input Parameter:
38: path - The dynamic library path, or PETSC_NULL
40: Level: developer
42: .seealso: SlepcInitialize()
43: @*/
44: PetscErrorCode EPSInitializePackage(char *path) {
45: static PetscTruth initialized = PETSC_FALSE;
46: char logList[256];
47: char *className;
48: PetscTruth opt;
52: if (initialized) return(0);
53: initialized = PETSC_TRUE;
54: /* Register Classes */
55: PetscCookieRegister("Eigenproblem Solver",&EPS_COOKIE);
56: /* Register Constructors */
57: EPSRegisterAll(path);
58: /* Register Events */
59: PetscLogEventRegister("EPSSetUp",EPS_COOKIE,&EPS_SetUp);
60: PetscLogEventRegister("EPSSolve",EPS_COOKIE,&EPS_Solve);
61: PetscLogEventRegister("EPSDense",EPS_COOKIE,&EPS_Dense);
62: /* Process info exclusions */
63: PetscOptionsGetString(PETSC_NULL, "-log_info_exclude", logList, 256, &opt);
64: if (opt) {
65: PetscStrstr(logList, "eps", &className);
66: if (className) {
67: PetscInfoDeactivateClass(EPS_COOKIE);
68: }
69: }
70: /* Process summary exclusions */
71: PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);
72: if (opt) {
73: PetscStrstr(logList, "eps", &className);
74: if (className) {
75: PetscLogEventDeactivateClass(EPS_COOKIE);
76: }
77: }
78: return(0);
79: }
83: /*@C
84: EPSView - Prints the EPS data structure.
86: Collective on EPS
88: Input Parameters:
89: + eps - the eigenproblem solver context
90: - viewer - optional visualization context
92: Options Database Key:
93: . -eps_view - Calls EPSView() at end of EPSSolve()
95: Note:
96: The available visualization contexts include
97: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
98: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
99: output where only the first processor opens
100: the file. All other processors send their
101: data to the first processor to print.
103: The user can open an alternative visualization context with
104: PetscViewerASCIIOpen() - output to a specified file.
106: Level: beginner
108: .seealso: STView(), PetscViewerASCIIOpen()
109: @*/
110: PetscErrorCode EPSView(EPS eps,PetscViewer viewer)
111: {
113: const EPSType type;
114: const char *extr, *which;
115: PetscTruth isascii;
119: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(((PetscObject)eps)->comm);
123: #if defined(PETSC_USE_COMPLEX)
124: #define HERM "hermitian"
125: #else
126: #define HERM "symmetric"
127: #endif
128: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
129: if (isascii) {
130: PetscViewerASCIIPrintf(viewer,"EPS Object:\n");
131: switch (eps->problem_type) {
132: case EPS_HEP: type = HERM " eigenvalue problem"; break;
133: case EPS_GHEP: type = "generalized " HERM " eigenvalue problem"; break;
134: case EPS_NHEP: type = "non-" HERM " eigenvalue problem"; break;
135: case EPS_GNHEP: type = "generalized non-" HERM " eigenvalue problem"; break;
136: case EPS_PGNHEP: type = "generalized non-" HERM " eigenvalue problem with " HERM " positive definite B"; break;
137: case 0: type = "not yet set"; break;
138: default: SETERRQ(1,"Wrong value of eps->problem_type");
139: }
140: PetscViewerASCIIPrintf(viewer," problem type: %s\n",type);
141: EPSGetType(eps,&type);
142: if (type) {
143: PetscViewerASCIIPrintf(viewer," method: %s",type);
144: switch (eps->solverclass) {
145: case EPS_ONE_SIDE:
146: PetscViewerASCIIPrintf(viewer,"\n",type); break;
147: case EPS_TWO_SIDE:
148: PetscViewerASCIIPrintf(viewer," (two-sided)\n",type); break;
149: default: SETERRQ(1,"Wrong value of eps->solverclass");
150: }
151: } else {
152: PetscViewerASCIIPrintf(viewer," method: not yet set\n");
153: }
154: if (eps->ops->view) {
155: PetscViewerASCIIPushTab(viewer);
156: (*eps->ops->view)(eps,viewer);
157: PetscViewerASCIIPopTab(viewer);
158: }
159: if (eps->extraction) {
160: switch (eps->extraction) {
161: case EPS_RITZ: extr = "Rayleigh-Ritz"; break;
162: case EPS_HARMONIC: extr = "harmonic Ritz"; break;
163: case EPS_REFINED: extr = "refined Ritz"; break;
164: case EPS_REFINED_HARMONIC: extr = "refined harmonic Ritz"; break;
165: default: SETERRQ(1,"Wrong value of eps->extraction");
166: }
167: PetscViewerASCIIPrintf(viewer," extraction type: %s\n",extr);
168: }
169: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: ");
170: if (eps->target_set) {
171: #if !defined(PETSC_USE_COMPLEX)
172: PetscViewerASCIIPrintf(viewer,"closest to target: %g\n",eps->target);
173: #else
174: PetscViewerASCIIPrintf(viewer,"closest to target: %g+%g i\n",PetscRealPart(eps->target),PetscImaginaryPart(eps->target));
175: #endif
176: } else {
177: switch (eps->which) {
178: case EPS_LARGEST_MAGNITUDE: which = "largest eigenvalues in magnitude"; break;
179: case EPS_SMALLEST_MAGNITUDE: which = "smallest eigenvalues in magnitude"; break;
180: case EPS_LARGEST_REAL: which = "largest real parts"; break;
181: case EPS_SMALLEST_REAL: which = "smallest real parts"; break;
182: case EPS_LARGEST_IMAGINARY: which = "largest imaginary parts"; break;
183: case EPS_SMALLEST_IMAGINARY: which = "smallest imaginary parts"; break;
184: default: SETERRQ(1,"Wrong value of eps->which");
185: }
186: PetscViewerASCIIPrintf(viewer,"%s\n",which);
187: }
188: PetscViewerASCIIPrintf(viewer," number of eigenvalues (nev): %d\n",eps->nev);
189: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %d\n",eps->ncv);
190: PetscViewerASCIIPrintf(viewer," maximum dimension of projected problem (mpd): %d\n",eps->mpd);
191: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %d\n", eps->max_it);
192: PetscViewerASCIIPrintf(viewer," tolerance: %g\n",eps->tol);
193: PetscViewerASCIIPrintf(viewer," dimension of user-provided deflation space: %d\n",eps->nds);
194: PetscViewerASCIIPushTab(viewer);
195: IPView(eps->ip,viewer);
196: STView(eps->OP,viewer);
197: PetscViewerASCIIPopTab(viewer);
198: } else {
199: if (eps->ops->view) {
200: (*eps->ops->view)(eps,viewer);
201: }
202: STView(eps->OP,viewer);
203: }
204: return(0);
205: }
209: /*@C
210: EPSCreate - Creates the default EPS context.
212: Collective on MPI_Comm
214: Input Parameter:
215: . comm - MPI communicator
217: Output Parameter:
218: . eps - location to put the EPS context
220: Note:
221: The default EPS type is EPSKRYLOVSCHUR
223: Level: beginner
225: .seealso: EPSSetUp(), EPSSolve(), EPSDestroy(), EPS
226: @*/
227: PetscErrorCode EPSCreate(MPI_Comm comm,EPS *outeps)
228: {
230: EPS eps;
234: *outeps = 0;
236: PetscHeaderCreate(eps,_p_EPS,struct _EPSOps,EPS_COOKIE,-1,"EPS",comm,EPSDestroy,EPSView);
237: *outeps = eps;
239: PetscMemzero(eps->ops,sizeof(struct _EPSOps));
241: eps->max_it = 0;
242: eps->nev = 1;
243: eps->ncv = 0;
244: eps->mpd = 0;
245: eps->allocated_ncv = 0;
246: eps->nds = 0;
247: eps->tol = 1e-7;
248: eps->which = EPS_LARGEST_MAGNITUDE;
249: eps->target = 0.0;
250: eps->target_set = PETSC_FALSE;
251: eps->evecsavailable = PETSC_FALSE;
252: eps->problem_type = (EPSProblemType)0;
253: eps->extraction = (EPSExtraction)0;
254: eps->solverclass = (EPSClass)0;
256: eps->vec_initial = 0;
257: eps->vec_initial_left= 0;
258: eps->V = 0;
259: eps->AV = 0;
260: eps->W = 0;
261: eps->T = 0;
262: eps->DS = 0;
263: eps->ds_ortho = PETSC_TRUE;
264: eps->eigr = 0;
265: eps->eigi = 0;
266: eps->errest = 0;
267: eps->errest_left = 0;
268: eps->OP = 0;
269: eps->ip = 0;
270: eps->data = 0;
271: eps->nconv = 0;
272: eps->its = 0;
273: eps->perm = PETSC_NULL;
275: eps->nwork = 0;
276: eps->work = 0;
277: eps->isgeneralized = PETSC_FALSE;
278: eps->ishermitian = PETSC_FALSE;
279: eps->ispositive = PETSC_FALSE;
280: eps->setupcalled = 0;
281: eps->reason = EPS_CONVERGED_ITERATING;
283: eps->numbermonitors = 0;
285: STCreate(comm,&eps->OP);
286: PetscLogObjectParent(eps,eps->OP);
287: IPCreate(comm,&eps->ip);
288: IPSetOptionsPrefix(eps->ip,((PetscObject)eps)->prefix);
289: IPAppendOptionsPrefix(eps->ip,"eps_");
290: PetscLogObjectParent(eps,eps->ip);
291: PetscPublishAll(eps);
292: return(0);
293: }
294:
297: /*@C
298: EPSSetType - Selects the particular solver to be used in the EPS object.
300: Collective on EPS
302: Input Parameters:
303: + eps - the eigensolver context
304: - type - a known method
306: Options Database Key:
307: . -eps_type <method> - Sets the method; use -help for a list
308: of available methods
309:
310: Notes:
311: See "slepc/include/slepceps.h" for available methods. The default
312: is EPSKRYLOVSCHUR.
314: Normally, it is best to use the EPSSetFromOptions() command and
315: then set the EPS type from the options database rather than by using
316: this routine. Using the options database provides the user with
317: maximum flexibility in evaluating the different available methods.
318: The EPSSetType() routine is provided for those situations where it
319: is necessary to set the iterative solver independently of the command
320: line or options database.
322: Level: intermediate
324: .seealso: STSetType(), EPSType
325: @*/
326: PetscErrorCode EPSSetType(EPS eps,const EPSType type)
327: {
328: PetscErrorCode ierr,(*r)(EPS);
329: PetscTruth match;
335: PetscTypeCompare((PetscObject)eps,type,&match);
336: if (match) return(0);
338: if (eps->data) {
339: /* destroy the old private EPS context */
340: (*eps->ops->destroy)(eps);
341: eps->data = 0;
342: }
344: PetscFListFind(EPSList,((PetscObject)eps)->comm,type,(void (**)(void)) &r);
346: if (!r) SETERRQ1(1,"Unknown EPS type given: %s",type);
348: eps->setupcalled = 0;
349: PetscMemzero(eps->ops,sizeof(struct _EPSOps));
350: (*r)(eps);
352: PetscObjectChangeTypeName((PetscObject)eps,type);
353: return(0);
354: }
358: /*@C
359: EPSGetType - Gets the EPS type as a string from the EPS object.
361: Not Collective
363: Input Parameter:
364: . eps - the eigensolver context
366: Output Parameter:
367: . name - name of EPS method
369: Level: intermediate
371: .seealso: EPSSetType()
372: @*/
373: PetscErrorCode EPSGetType(EPS eps,const EPSType *type)
374: {
378: *type = ((PetscObject)eps)->type_name;
379: return(0);
380: }
382: /*MC
383: EPSRegisterDynamic - Adds a method to the eigenproblem solver package.
385: Synopsis:
386: EPSRegisterDynamic(char *name_solver,char *path,char *name_create,PetscErrorCode (*routine_create)(EPS))
388: Not Collective
390: Input Parameters:
391: + name_solver - name of a new user-defined solver
392: . path - path (either absolute or relative) the library containing this solver
393: . name_create - name of routine to create the solver context
394: - routine_create - routine to create the solver context
396: Notes:
397: EPSRegisterDynamic() may be called multiple times to add several user-defined solvers.
399: If dynamic libraries are used, then the fourth input argument (routine_create)
400: is ignored.
402: Sample usage:
403: .vb
404: EPSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
405: "MySolverCreate",MySolverCreate);
406: .ve
408: Then, your solver can be chosen with the procedural interface via
409: $ EPSSetType(eps,"my_solver")
410: or at runtime via the option
411: $ -eps_type my_solver
413: Level: advanced
415: Environmental variables such as ${PETSC_ARCH}, ${PETSC_DIR}, ${PETSC_LIB_DIR},
416: and others of the form ${any_environmental_variable} occuring in pathname will be
417: replaced with appropriate values.
419: .seealso: EPSRegisterDestroy(), EPSRegisterAll()
421: M*/
425: /*@C
426: EPSRegister - See EPSRegisterDynamic()
428: Level: advanced
429: @*/
430: PetscErrorCode EPSRegister(const char *sname,const char *path,const char *name,PetscErrorCode (*function)(EPS))
431: {
433: char fullname[256];
436: PetscFListConcat(path,name,fullname);
437: PetscFListAdd(&EPSList,sname,fullname,(void (*)(void))function);
438: return(0);
439: }
443: /*@
444: EPSRegisterDestroy - Frees the list of EPS methods that were
445: registered by EPSRegisterDynamic().
447: Not Collective
449: Level: advanced
451: .seealso: EPSRegisterDynamic(), EPSRegisterAll()
452: @*/
453: PetscErrorCode EPSRegisterDestroy(void)
454: {
458: PetscFListDestroy(&EPSList);
459: EPSRegisterAll(PETSC_NULL);
460: return(0);
461: }
465: /*@
466: EPSDestroy - Destroys the EPS context.
468: Collective on EPS
470: Input Parameter:
471: . eps - eigensolver context obtained from EPSCreate()
473: Level: beginner
475: .seealso: EPSCreate(), EPSSetUp(), EPSSolve()
476: @*/
477: PetscErrorCode EPSDestroy(EPS eps)
478: {
483: if (--((PetscObject)eps)->refct > 0) return(0);
485: /* if memory was published with AMS then destroy it */
486: PetscObjectDepublish(eps);
488: STDestroy(eps->OP);
489: IPDestroy(eps->ip);
491: if (eps->ops->destroy) {
492: (*eps->ops->destroy)(eps);
493: }
494:
495: PetscFree(eps->T);
496: PetscFree(eps->Tl);
497: PetscFree(eps->perm);
499: if (eps->vec_initial) {
500: VecDestroy(eps->vec_initial);
501: }
503: if (eps->vec_initial_left) {
504: VecDestroy(eps->vec_initial_left);
505: }
507: if (eps->nds > 0) {
508: VecDestroyVecs(eps->DS, eps->nds);
509: }
510:
511: PetscFree(eps->DSV);
513: EPSMonitorCancel(eps);
515: PetscHeaderDestroy(eps);
516: return(0);
517: }
521: /*@
522: EPSSetTarget - Sets the value of the target.
524: Not collective
526: Input Parameters:
527: + eps - eigensolver context
528: - target - the value of the target
530: Notes:
531: The target is a scalar value used to determine the portion of the spectrum
532: of interest.
534: If the target is not specified, then eigenvalues are computed according to
535: the which parameter (see EPSSetWhichEigenpairs()).
536:
537: If the target is specified, then the sought-after eigenvalues are those
538: closest to the target.
540: Level: beginner
542: .seealso: EPSGetTarget(), EPSSetWhichEigenpairs()
543: @*/
544: PetscErrorCode EPSSetTarget(EPS eps,PetscScalar target)
545: {
548: eps->target = target;
549: eps->target_set = PETSC_TRUE;
550: return(0);
551: }
555: /*@
556: EPSGetTarget - Gets the value of the target.
558: Not collective
560: Input Parameter:
561: . eps - eigensolver context
563: Output Parameter:
564: . target - the value of the target
566: Level: beginner
568: Note:
569: If the target was not set by the user, then zero is returned.
571: .seealso: EPSSetTarget()
572: @*/
573: PetscErrorCode EPSGetTarget(EPS eps,PetscScalar* target)
574: {
577: if (target) {
578: if (eps->target_set) *target = eps->target;
579: else *target = 0.0;
580: }
581: return(0);
582: }
586: /*@
587: EPSSetST - Associates a spectral transformation object to the
588: eigensolver.
590: Collective on EPS
592: Input Parameters:
593: + eps - eigensolver context obtained from EPSCreate()
594: - st - the spectral transformation object
596: Note:
597: Use EPSGetST() to retrieve the spectral transformation context (for example,
598: to free it at the end of the computations).
600: Level: advanced
602: .seealso: EPSGetST()
603: @*/
604: PetscErrorCode EPSSetST(EPS eps,ST st)
605: {
612: PetscObjectReference((PetscObject)st);
613: STDestroy(eps->OP);
614: eps->OP = st;
615: return(0);
616: }
620: /*@C
621: EPSGetST - Obtain the spectral transformation (ST) object associated
622: to the eigensolver object.
624: Not Collective
626: Input Parameters:
627: . eps - eigensolver context obtained from EPSCreate()
629: Output Parameter:
630: . st - spectral transformation context
632: Level: beginner
634: .seealso: EPSSetST()
635: @*/
636: PetscErrorCode EPSGetST(EPS eps, ST *st)
637: {
641: *st = eps->OP;
642: return(0);
643: }
647: /*@
648: EPSSetIP - Associates an inner product object to the
649: eigensolver.
651: Collective on EPS
653: Input Parameters:
654: + eps - eigensolver context obtained from EPSCreate()
655: - ip - the inner product object
657: Note:
658: Use EPSGetIP() to retrieve the inner product context (for example,
659: to free it at the end of the computations).
661: Level: advanced
663: .seealso: EPSGetIP()
664: @*/
665: PetscErrorCode EPSSetIP(EPS eps,IP ip)
666: {
673: PetscObjectReference((PetscObject)ip);
674: IPDestroy(eps->ip);
675: eps->ip = ip;
676: return(0);
677: }
681: /*@C
682: EPSGetIP - Obtain the inner product object associated
683: to the eigensolver object.
685: Not Collective
687: Input Parameters:
688: . eps - eigensolver context obtained from EPSCreate()
690: Output Parameter:
691: . ip - inner product context
693: Level: advanced
695: .seealso: EPSSetIP()
696: @*/
697: PetscErrorCode EPSGetIP(EPS eps,IP *ip)
698: {
702: *ip = eps->ip;
703: return(0);
704: }
708: /*@
709: EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
710: eigenvalue problem.
712: Not collective
714: Input Parameter:
715: . eps - the eigenproblem solver context
717: Output Parameter:
718: . is - the answer
720: Level: intermediate
722: @*/
723: PetscErrorCode EPSIsGeneralized(EPS eps,PetscTruth* is)
724: {
726: Mat B;
730: STGetOperators(eps->OP,PETSC_NULL,&B);
731: if( B ) *is = PETSC_TRUE;
732: else *is = PETSC_FALSE;
733: if( eps->setupcalled ) {
734: if( eps->isgeneralized != *is ) {
735: SETERRQ(0,"Warning: Inconsistent EPS state");
736: }
737: }
738: return(0);
739: }
743: /*@
744: EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
745: eigenvalue problem.
747: Not collective
749: Input Parameter:
750: . eps - the eigenproblem solver context
752: Output Parameter:
753: . is - the answer
755: Level: intermediate
757: @*/
758: PetscErrorCode EPSIsHermitian(EPS eps,PetscTruth* is)
759: {
762: if( eps->ishermitian ) *is = PETSC_TRUE;
763: else *is = PETSC_FALSE;
764: return(0);
765: }