Actual source code: basic.c
1: /*
2: The basic EPS routines, Create, View, etc. are here.
3: */
4: #include src/eps/epsimpl.h
6: PetscFList EPSList = 0;
7: PetscCookie EPS_COOKIE = 0;
8: PetscEvent EPS_SetUp = 0, EPS_Solve = 0, EPS_Orthogonalize = 0;
12: /*@C
13: EPSInitializePackage - This function initializes everything in the EPS package. It is called
14: from PetscDLLibraryRegister() when using dynamic libraries, and on the first call to EPSCreate()
15: when using static libraries.
17: Input Parameter:
18: path - The dynamic library path, or PETSC_NULL
20: Level: developer
22: .seealso: SlepcInitialize()
23: @*/
24: PetscErrorCode EPSInitializePackage(char *path) {
25: static PetscTruth initialized = PETSC_FALSE;
26: char logList[256];
27: char *className;
28: PetscTruth opt;
32: if (initialized) return(0);
33: initialized = PETSC_TRUE;
34: /* Register Classes */
35: PetscLogClassRegister(&EPS_COOKIE,"Eigenproblem Solver");
36: /* Register Constructors */
37: EPSRegisterAll(path);
38: /* Register Events */
39: PetscLogEventRegister(&EPS_SetUp,"EPSSetUp",EPS_COOKIE);
40: PetscLogEventRegister(&EPS_Solve,"EPSSolve",EPS_COOKIE);
41: PetscLogEventRegister(&EPS_Orthogonalize,"EPSOrthogonalize",EPS_COOKIE);
42: /* Process info exclusions */
43: PetscOptionsGetString(PETSC_NULL, "-log_info_exclude", logList, 256, &opt);
44: if (opt) {
45: PetscStrstr(logList, "eps", &className);
46: if (className) {
47: PetscInfoDeactivateClass(EPS_COOKIE);
48: }
49: }
50: /* Process summary exclusions */
51: PetscOptionsGetString(PETSC_NULL, "-log_summary_exclude", logList, 256, &opt);
52: if (opt) {
53: PetscStrstr(logList, "eps", &className);
54: if (className) {
55: PetscLogEventDeactivateClass(EPS_COOKIE);
56: }
57: }
58: return(0);
59: }
63: /*@C
64: EPSView - Prints the EPS data structure.
66: Collective on EPS
68: Input Parameters:
69: + eps - the eigenproblem solver context
70: - viewer - optional visualization context
72: Options Database Key:
73: . -eps_view - Calls EPSView() at end of EPSSolve()
75: Note:
76: The available visualization contexts include
77: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
78: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
79: output where only the first processor opens
80: the file. All other processors send their
81: data to the first processor to print.
83: The user can open an alternative visualization context with
84: PetscViewerASCIIOpen() - output to a specified file.
86: Level: beginner
88: .seealso: STView(), PetscViewerASCIIOpen()
89: @*/
90: PetscErrorCode EPSView(EPS eps,PetscViewer viewer)
91: {
93: const char *type, *which;
94: PetscTruth isascii;
98: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(eps->comm);
102: #if defined(PETSC_USE_COMPLEX)
103: #define HERM "hermitian"
104: #else
105: #define HERM "symmetric"
106: #endif
107: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
108: if (isascii) {
109: PetscViewerASCIIPrintf(viewer,"EPS Object:\n");
110: switch (eps->problem_type) {
111: case EPS_HEP: type = HERM " eigenvalue problem"; break;
112: case EPS_GHEP: type = "generalized " HERM " eigenvalue problem"; break;
113: case EPS_NHEP: type = "non-" HERM " eigenvalue problem"; break;
114: case EPS_GNHEP: type = "generalized non-" HERM " eigenvalue problem"; break;
115: case 0: type = "not yet set"; break;
116: default: SETERRQ(1,"Wrong value of eps->problem_type");
117: }
118: PetscViewerASCIIPrintf(viewer," problem type: %s\n",type);
119: EPSGetType(eps,&type);
120: if (type) {
121: PetscViewerASCIIPrintf(viewer," method: %s",type);
122: switch (eps->solverclass) {
123: case EPS_ONE_SIDE:
124: PetscViewerASCIIPrintf(viewer,"\n",type); break;
125: case EPS_TWO_SIDE:
126: PetscViewerASCIIPrintf(viewer," (two-sided)\n",type); break;
127: default: SETERRQ(1,"Wrong value of eps->solverclass");
128: }
129: } else {
130: PetscViewerASCIIPrintf(viewer," method: not yet set\n");
131: }
132: if (eps->ops->view) {
133: PetscViewerASCIIPushTab(viewer);
134: (*eps->ops->view)(eps,viewer);
135: PetscViewerASCIIPopTab(viewer);
136: }
137: switch (eps->which) {
138: case EPS_LARGEST_MAGNITUDE: which = "largest eigenvalues in magnitude"; break;
139: case EPS_SMALLEST_MAGNITUDE: which = "smallest eigenvalues in magnitude"; break;
140: case EPS_LARGEST_REAL: which = "largest real parts"; break;
141: case EPS_SMALLEST_REAL: which = "smallest real parts"; break;
142: case EPS_LARGEST_IMAGINARY: which = "largest imaginary parts"; break;
143: case EPS_SMALLEST_IMAGINARY: which = "smallest imaginary parts"; break;
144: default: SETERRQ(1,"Wrong value of eps->which");
145: }
146: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: %s\n",which);
147: PetscViewerASCIIPrintf(viewer," number of eigenvalues (nev): %d\n",eps->nev);
148: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %d\n",eps->ncv);
149: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %d\n", eps->max_it);
150: PetscViewerASCIIPrintf(viewer," tolerance: %g\n",eps->tol);
151: PetscViewerASCIIPrintf(viewer," orthogonalization method: ");
152: switch (eps->orthog_type) {
153: case EPS_MGS_ORTH:
154: PetscViewerASCIIPrintf(viewer,"modified Gram-Schmidt\n");
155: break;
156: case EPS_CGS_ORTH:
157: PetscViewerASCIIPrintf(viewer,"classical Gram-Schmidt\n");
158: break;
159: default: SETERRQ(1,"Wrong value of eps->orth_type");
160: }
161: PetscViewerASCIIPrintf(viewer," orthogonalization refinement: ");
162: switch (eps->orthog_ref) {
163: case EPS_ORTH_REFINE_NEVER:
164: PetscViewerASCIIPrintf(viewer,"never\n");
165: break;
166: case EPS_ORTH_REFINE_IFNEEDED:
167: PetscViewerASCIIPrintf(viewer,"if needed (eta: %f)\n",eps->orthog_eta);
168: break;
169: case EPS_ORTH_REFINE_ALWAYS:
170: PetscViewerASCIIPrintf(viewer,"always\n");
171: break;
172: default: SETERRQ(1,"Wrong value of eps->orth_ref");
173: }
174: PetscViewerASCIIPrintf(viewer," dimension of user-provided deflation space: %d\n",eps->nds);
175: PetscViewerASCIIPushTab(viewer);
176: STView(eps->OP,viewer);
177: PetscViewerASCIIPopTab(viewer);
178: } else {
179: if (eps->ops->view) {
180: (*eps->ops->view)(eps,viewer);
181: }
182: STView(eps->OP,viewer);
183: }
184: return(0);
185: }
189: static PetscErrorCode EPSPublish_Petsc(PetscObject object)
190: {
192: return(0);
193: }
197: /*@C
198: EPSCreate - Creates the default EPS context.
200: Collective on MPI_Comm
202: Input Parameter:
203: . comm - MPI communicator
205: Output Parameter:
206: . eps - location to put the EPS context
208: Note:
209: The default EPS type is EPSKRYLOVSCHUR
211: Level: beginner
213: .seealso: EPSSetUp(), EPSSolve(), EPSDestroy(), EPS
214: @*/
215: PetscErrorCode EPSCreate(MPI_Comm comm,EPS *outeps)
216: {
218: EPS eps;
222: *outeps = 0;
224: PetscHeaderCreate(eps,_p_EPS,struct _EPSOps,EPS_COOKIE,-1,"EPS",comm,EPSDestroy,EPSView);
225: PetscLogObjectCreate(eps);
226: *outeps = eps;
228: eps->bops->publish = EPSPublish_Petsc;
229: PetscMemzero(eps->ops,sizeof(struct _EPSOps));
231: eps->type = -1;
232: eps->max_it = 0;
233: eps->nev = 1;
234: eps->ncv = 0;
235: eps->allocated_ncv = 0;
236: eps->nds = 0;
237: eps->tol = 0.0;
238: eps->which = EPS_LARGEST_MAGNITUDE;
239: eps->evecsavailable = PETSC_FALSE;
240: eps->problem_type = (EPSProblemType)0;
241: eps->solverclass = (EPSClass)0;
243: eps->vec_initial = 0;
244: eps->vec_initial_left= 0;
245: eps->V = 0;
246: eps->AV = 0;
247: eps->W = 0;
248: eps->AW = 0;
249: eps->T = 0;
250: eps->DS = 0;
251: eps->ds_ortho = PETSC_TRUE;
252: eps->eigr = 0;
253: eps->eigi = 0;
254: eps->errest = 0;
255: eps->errest_left = 0;
256: eps->OP = 0;
257: eps->data = 0;
258: eps->nconv = 0;
259: eps->its = 0;
260: eps->perm = PETSC_NULL;
262: eps->nwork = 0;
263: eps->work = 0;
264: eps->isgeneralized = PETSC_FALSE;
265: eps->ishermitian = PETSC_FALSE;
266: eps->setupcalled = 0;
267: eps->reason = EPS_CONVERGED_ITERATING;
269: eps->numbermonitors = 0;
271: eps->orthog_type = EPS_CGS_ORTH;
272: eps->orthog_ref = EPS_ORTH_REFINE_IFNEEDED;
273: eps->orthog_eta = PETSC_DEFAULT;
275: STCreate(comm,&eps->OP);
276: PetscLogObjectParent(eps,eps->OP);
277: PetscPublishAll(eps);
278: return(0);
279: }
280:
283: /*@C
284: EPSSetType - Selects the particular solver to be used in the EPS object.
286: Collective on EPS
288: Input Parameters:
289: + eps - the eigensolver context
290: - type - a known method
292: Options Database Key:
293: . -eps_type <method> - Sets the method; use -help for a list
294: of available methods
295:
296: Notes:
297: See "slepc/include/slepceps.h" for available methods. The default
298: is EPSKRYLOVSCHUR.
300: Normally, it is best to use the EPSSetFromOptions() command and
301: then set the EPS type from the options database rather than by using
302: this routine. Using the options database provides the user with
303: maximum flexibility in evaluating the different available methods.
304: The EPSSetType() routine is provided for those situations where it
305: is necessary to set the iterative solver independently of the command
306: line or options database.
308: Level: intermediate
310: .seealso: STSetType(), EPSType
311: @*/
312: PetscErrorCode EPSSetType(EPS eps,EPSType type)
313: {
314: PetscErrorCode ierr,(*r)(EPS);
315: PetscTruth match;
321: PetscTypeCompare((PetscObject)eps,type,&match);
322: if (match) return(0);
324: if (eps->data) {
325: /* destroy the old private EPS context */
326: (*eps->ops->destroy)(eps);
327: eps->data = 0;
328: }
330: PetscFListFind(eps->comm,EPSList,type,(void (**)(void)) &r);
332: if (!r) SETERRQ1(1,"Unknown EPS type given: %s",type);
334: eps->setupcalled = 0;
335: PetscMemzero(eps->ops,sizeof(struct _EPSOps));
336: (*r)(eps);
338: PetscObjectChangeTypeName((PetscObject)eps,type);
339: return(0);
340: }
344: /*@C
345: EPSGetType - Gets the EPS type as a string from the EPS object.
347: Not Collective
349: Input Parameter:
350: . eps - the eigensolver context
352: Output Parameter:
353: . name - name of EPS method
355: Level: intermediate
357: .seealso: EPSSetType()
358: @*/
359: PetscErrorCode EPSGetType(EPS eps,EPSType *type)
360: {
363: *type = eps->type_name;
364: return(0);
365: }
367: /*MC
368: EPSRegisterDynamic - Adds a method to the eigenproblem solver package.
370: Synopsis:
371: EPSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(EPS))
373: Not Collective
375: Input Parameters:
376: + name_solver - name of a new user-defined solver
377: . path - path (either absolute or relative) the library containing this solver
378: . name_create - name of routine to create the solver context
379: - routine_create - routine to create the solver context
381: Notes:
382: EPSRegisterDynamic() may be called multiple times to add several user-defined solvers.
384: If dynamic libraries are used, then the fourth input argument (routine_create)
385: is ignored.
387: Sample usage:
388: .vb
389: EPSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
390: "MySolverCreate",MySolverCreate);
391: .ve
393: Then, your solver can be chosen with the procedural interface via
394: $ EPSSetType(eps,"my_solver")
395: or at runtime via the option
396: $ -eps_type my_solver
398: Level: advanced
400: Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR}, ${BOPT},
401: and others of the form ${any_environmental_variable} occuring in pathname will be
402: replaced with appropriate values.
404: .seealso: EPSRegisterAll()
406: M*/
410: PetscErrorCode EPSRegister(const char *sname,const char *path,const char *name,int (*function)(EPS))
411: {
413: char fullname[256];
416: PetscFListConcat(path,name,fullname);
417: PetscFListAdd(&EPSList,sname,fullname,(void (*)(void))function);
418: return(0);
419: }
423: /*@
424: EPSDestroy - Destroys the EPS context.
426: Collective on EPS
428: Input Parameter:
429: . eps - eigensolver context obtained from EPSCreate()
431: Level: beginner
433: .seealso: EPSCreate(), EPSSetUp(), EPSSolve()
434: @*/
435: PetscErrorCode EPSDestroy(EPS eps)
436: {
441: if (--eps->refct > 0) return(0);
443: /* if memory was published with AMS then destroy it */
444: PetscObjectDepublish(eps);
446: STDestroy(eps->OP);
448: if (eps->ops->destroy) {
449: (*eps->ops->destroy)(eps);
450: }
451:
452: PetscFree(eps->T);
453: PetscFree(eps->Tl);
454: PetscFree(eps->perm);
456: if (eps->vec_initial) {
457: VecDestroy(eps->vec_initial);
458: }
460: if (eps->vec_initial_left) {
461: VecDestroy(eps->vec_initial_left);
462: }
464: if (eps->nds > 0) {
465: VecDestroyVecs(eps->DS, eps->nds);
466: }
467:
468: PetscFree(eps->DSV);
470: PetscLogObjectDestroy(eps);
471: PetscHeaderDestroy(eps);
472: return(0);
473: }
477: /*@
478: EPSSetST - Associates a spectral transformation object to the
479: eigensolver.
481: Collective on EPS
483: Input Parameters:
484: + eps - eigensolver context obtained from EPSCreate()
485: - st - the spectral transformation object
487: Note:
488: Use EPSGetST() to retrieve the spectral transformation context (for example,
489: to free it at the end of the computations).
491: Level: advanced
493: .seealso: EPSGetST()
494: @*/
495: PetscErrorCode EPSSetST(EPS eps,ST st)
496: {
503: STDestroy(eps->OP);
504: eps->OP = st;
505: PetscObjectReference((PetscObject)eps->OP);
506: return(0);
507: }
511: /*@C
512: EPSGetST - Obtain the spectral transformation (ST) object associated
513: to the eigensolver object.
515: Not Collective
517: Input Parameters:
518: . eps - eigensolver context obtained from EPSCreate()
520: Output Parameter:
521: . st - spectral transformation context
523: Level: beginner
525: .seealso: EPSSetST()
526: @*/
527: PetscErrorCode EPSGetST(EPS eps, ST *st)
528: {
531: *st = eps->OP;
532: return(0);
533: }
537: /*@
538: EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
539: eigenvalue problem.
541: Not collective
543: Input Parameter:
544: . eps - the eigenproblem solver context
546: Output Parameter:
547: . is - the answer
549: Level: intermediate
551: @*/
552: PetscErrorCode EPSIsGeneralized(EPS eps,PetscTruth* is)
553: {
555: Mat B;
559: STGetOperators(eps->OP,PETSC_NULL,&B);
560: if( B ) *is = PETSC_TRUE;
561: else *is = PETSC_FALSE;
562: if( eps->setupcalled ) {
563: if( eps->isgeneralized != *is ) {
564: SETERRQ(0,"Warning: Inconsistent EPS state");
565: }
566: }
567: return(0);
568: }
572: /*@
573: EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
574: eigenvalue problem.
576: Not collective
578: Input Parameter:
579: . eps - the eigenproblem solver context
581: Output Parameter:
582: . is - the answer
584: Level: intermediate
586: @*/
587: PetscErrorCode EPSIsHermitian(EPS eps,PetscTruth* is)
588: {
591: if( eps->ishermitian ) *is = PETSC_TRUE;
592: else *is = PETSC_FALSE;
593: return(0);
594: }