Actual source code: basic.c
1: /*
2: The basic EPS routines, Create, View, etc. are here.
3: */
4: #include src/eps/epsimpl.h
6: PetscTruth EPSRegisterAllCalled = PETSC_FALSE;
7: PetscFList EPSList = 0;
11: /*@C
12: EPSView - Prints the EPS data structure.
14: Collective on EPS
16: Input Parameters:
17: + eps - the eigenproblem solver context
18: - viewer - optional visualization context
20: Options Database Key:
21: . -eps_view - Calls EPSView() at end of EPSSolve()
23: Note:
24: The available visualization contexts include
25: + PETSC_VIEWER_STDOUT_SELF - standard output (default)
26: - PETSC_VIEWER_STDOUT_WORLD - synchronized standard
27: output where only the first processor opens
28: the file. All other processors send their
29: data to the first processor to print.
31: The user can open an alternative visualization context with
32: PetscViewerASCIIOpen() - output to a specified file.
34: Level: beginner
36: .seealso: STView(), PetscViewerASCIIOpen()
37: @*/
38: PetscErrorCode EPSView(EPS eps,PetscViewer viewer)
39: {
40: PetscErrorCode ierr;
41: char *type, *which;
42: PetscTruth isascii;
46: if (!viewer) viewer = PETSC_VIEWER_STDOUT_(eps->comm);
50: #if defined(PETSC_USE_COMPLEX)
51: #define HERM "hermitian"
52: #else
53: #define HERM "symmetric"
54: #endif
55: PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&isascii);
56: if (isascii) {
57: PetscViewerASCIIPrintf(viewer,"EPS Object:\n");
58: switch (eps->problem_type) {
59: case EPS_HEP: type = HERM " eigenvalue problem"; break;
60: case EPS_GHEP: type = "generalized " HERM " eigenvalue problem"; break;
61: case EPS_NHEP: type = "non-" HERM " eigenvalue problem"; break;
62: case EPS_GNHEP: type = "generalized non-" HERM " eigenvalue problem"; break;
63: case 0: type = "not yet set"; break;
64: default: SETERRQ(1,"Wrong value of eps->problem_type");
65: }
66: PetscViewerASCIIPrintf(viewer," problem type: %s\n",type);
67: EPSGetType(eps,&type);
68: if (type) {
69: PetscViewerASCIIPrintf(viewer," method: %s\n",type);
70: } else {
71: PetscViewerASCIIPrintf(viewer," method: not yet set\n");
72: }
73: if (eps->ops->view) {
74: PetscViewerASCIIPushTab(viewer);
75: (*eps->ops->view)(eps,viewer);
76: PetscViewerASCIIPopTab(viewer);
77: }
78: switch (eps->which) {
79: case EPS_LARGEST_MAGNITUDE: which = "largest eigenvalues in magnitude"; break;
80: case EPS_SMALLEST_MAGNITUDE: which = "smallest eigenvalues in magnitude"; break;
81: case EPS_LARGEST_REAL: which = "largest real parts"; break;
82: case EPS_SMALLEST_REAL: which = "smallest real parts"; break;
83: case EPS_LARGEST_IMAGINARY: which = "largest imaginary parts"; break;
84: case EPS_SMALLEST_IMAGINARY: which = "smallest imaginary parts"; break;
85: default: SETERRQ(1,"Wrong value of eps->which");
86: }
87: PetscViewerASCIIPrintf(viewer," selected portion of the spectrum: %s\n",which);
88: PetscViewerASCIIPrintf(viewer," number of eigenvalues (nev): %d\n",eps->nev);
89: PetscViewerASCIIPrintf(viewer," number of column vectors (ncv): %d\n",eps->ncv);
90: PetscViewerASCIIPrintf(viewer," maximum number of iterations: %d\n", eps->max_it);
91: PetscViewerASCIIPrintf(viewer," tolerance: %g\n",eps->tol);
92: PetscViewerASCIIPrintf(viewer," orthogonalization method: ");
93: switch (eps->orthog_type) {
94: case EPS_MGS_ORTH:
95: PetscViewerASCIIPrintf(viewer,"modified Gram-Schmidt\n");
96: break;
97: case EPS_CGS_ORTH:
98: PetscViewerASCIIPrintf(viewer,"classical Gram-Schmidt\n");
99: break;
100: default: SETERRQ(1,"Wrong value of eps->orth_type");
101: }
102: PetscViewerASCIIPrintf(viewer," orthogonalization refinement: ");
103: switch (eps->orthog_ref) {
104: case EPS_ORTH_REFINE_NEVER:
105: PetscViewerASCIIPrintf(viewer,"never\n");
106: break;
107: case EPS_ORTH_REFINE_IFNEEDED:
108: PetscViewerASCIIPrintf(viewer,"if needed (eta: %f)\n",eps->orthog_eta);
109: break;
110: case EPS_ORTH_REFINE_ALWAYS:
111: PetscViewerASCIIPrintf(viewer,"always\n");
112: break;
113: default: SETERRQ(1,"Wrong value of eps->orth_ref");
114: }
115: PetscViewerASCIIPrintf(viewer," dimension of user-provided deflation space: %d\n",eps->nds);
116: PetscViewerASCIIPushTab(viewer);
117: STView(eps->OP,viewer);
118: PetscViewerASCIIPopTab(viewer);
119: } else {
120: if (eps->ops->view) {
121: (*eps->ops->view)(eps,viewer);
122: }
123: STView(eps->OP,viewer);
124: }
125: return(0);
126: }
130: static PetscErrorCode EPSPublish_Petsc(PetscObject object)
131: {
133: return(0);
134: }
138: /*@C
139: EPSCreate - Creates the default EPS context.
141: Collective on MPI_Comm
143: Input Parameter:
144: . comm - MPI communicator
146: Output Parameter:
147: . eps - location to put the EPS context
149: Note:
150: The default EPS type is EPSARNOLDI
152: Level: beginner
154: .seealso: EPSSetUp(), EPSSolve(), EPSDestroy(), EPS
155: @*/
156: PetscErrorCode EPSCreate(MPI_Comm comm,EPS *outeps)
157: {
158: PetscErrorCode ierr;
159: EPS eps;
163: *outeps = 0;
164: #ifndef PETSC_USE_DYNAMIC_LIBRARIES
165: EPSRegisterAll(PETSC_NULL);
166: #endif
168: PetscHeaderCreate(eps,_p_EPS,struct _EPSOps,EPS_COOKIE,-1,"EPS",comm,EPSDestroy,EPSView);
169: PetscLogObjectCreate(eps);
170: *outeps = eps;
172: eps->bops->publish = EPSPublish_Petsc;
173: PetscMemzero(eps->ops,sizeof(struct _EPSOps));
175: eps->type = -1;
176: eps->max_it = 0;
177: eps->nev = 1;
178: eps->ncv = 0;
179: eps->allocated_ncv = 0;
180: eps->nds = 0;
181: eps->tol = 0.0;
182: eps->which = EPS_LARGEST_MAGNITUDE;
183: eps->evecsavailable = PETSC_FALSE;
184: eps->problem_type = (EPSProblemType)0;
186: eps->vec_initial = 0;
187: eps->vec_initial_set = PETSC_FALSE;
188: eps->V = 0;
189: eps->AV = 0;
190: eps->T = 0;
191: eps->DS = 0;
192: eps->ds_ortho = PETSC_TRUE;
193: eps->eigr = 0;
194: eps->eigi = 0;
195: eps->errest = 0;
196: eps->OP = 0;
197: eps->data = 0;
198: eps->nconv = 0;
199: eps->its = 0;
200: eps->perm = PETSC_NULL;
202: eps->nwork = 0;
203: eps->work = 0;
204: eps->isgeneralized = PETSC_FALSE;
205: eps->ishermitian = PETSC_FALSE;
206: eps->setupcalled = 0;
207: eps->reason = EPS_CONVERGED_ITERATING;
209: eps->numbermonitors = 0;
211: eps->orthog_type = EPS_CGS_ORTH;
212: eps->orthog_ref = EPS_ORTH_REFINE_IFNEEDED;
213: eps->orthog_eta = 0.7071;
215: STCreate(comm,&eps->OP);
216: PetscLogObjectParent(eps,eps->OP);
217: PetscPublishAll(eps);
218: return(0);
219: }
220:
223: /*@C
224: EPSSetType - Selects the particular solver to be used in the EPS object.
226: Collective on EPS
228: Input Parameters:
229: + eps - the eigensolver context
230: - type - a known method
232: Options Database Key:
233: . -eps_type <method> - Sets the method; use -help for a list
234: of available methods
235:
236: Notes:
237: See "slepc/include/slepceps.h" for available methods. The default
238: is EPSARNOLDI.
240: Normally, it is best to use the EPSSetFromOptions() command and
241: then set the EPS type from the options database rather than by using
242: this routine. Using the options database provides the user with
243: maximum flexibility in evaluating the different available methods.
244: The EPSSetType() routine is provided for those situations where it
245: is necessary to set the iterative solver independently of the command
246: line or options database.
248: Level: intermediate
250: .seealso: STSetType(), EPSType
251: @*/
252: PetscErrorCode EPSSetType(EPS eps,EPSType type)
253: {
254: PetscErrorCode ierr,(*r)(EPS);
255: PetscTruth match;
261: PetscTypeCompare((PetscObject)eps,type,&match);
262: if (match) return(0);
264: if (eps->data) {
265: /* destroy the old private EPS context */
266: (*eps->ops->destroy)(eps);
267: eps->data = 0;
268: }
269: /* Get the function pointers for the iterative method requested */
270: if (!EPSRegisterAllCalled) {EPSRegisterAll(PETSC_NULL); }
272: PetscFListFind(eps->comm,EPSList,type,(void (**)(void)) &r);
274: if (!r) SETERRQ1(1,"Unknown EPS type given: %s",type);
276: eps->setupcalled = 0;
277: PetscMemzero(eps->ops,sizeof(struct _EPSOps));
278: (*r)(eps);
280: PetscObjectChangeTypeName((PetscObject)eps,type);
281: return(0);
282: }
286: /*@C
287: EPSRegisterDestroy - Frees the list of EPS methods that were
288: registered by EPSRegisterDynamic().
290: Not Collective
292: Level: advanced
294: .seealso: EPSRegisterDynamic(), EPSRegisterAll()
295: @*/
296: PetscErrorCode EPSRegisterDestroy(void)
297: {
298: PetscErrorCode ierr;
301: if (EPSList) {
302: PetscFListDestroy(&EPSList);
303: EPSList = 0;
304: }
305: EPSRegisterAllCalled = PETSC_FALSE;
306: return(0);
307: }
311: /*@C
312: EPSGetType - Gets the EPS type as a string from the EPS object.
314: Not Collective
316: Input Parameter:
317: . eps - the eigensolver context
319: Output Parameter:
320: . name - name of EPS method
322: Level: intermediate
324: .seealso: EPSSetType()
325: @*/
326: PetscErrorCode EPSGetType(EPS eps,EPSType *type)
327: {
330: *type = eps->type_name;
331: return(0);
332: }
334: /*MC
335: EPSRegisterDynamic - Adds a method to the eigenproblem solver package.
337: Synopsis:
338: EPSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(EPS))
340: Not Collective
342: Input Parameters:
343: + name_solver - name of a new user-defined solver
344: . path - path (either absolute or relative) the library containing this solver
345: . name_create - name of routine to create the solver context
346: - routine_create - routine to create the solver context
348: Notes:
349: EPSRegisterDynamic() may be called multiple times to add several user-defined solvers.
351: If dynamic libraries are used, then the fourth input argument (routine_create)
352: is ignored.
354: Sample usage:
355: .vb
356: EPSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
357: "MySolverCreate",MySolverCreate);
358: .ve
360: Then, your solver can be chosen with the procedural interface via
361: $ EPSSetType(eps,"my_solver")
362: or at runtime via the option
363: $ -eps_type my_solver
365: Level: advanced
367: Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR}, ${BOPT},
368: and others of the form ${any_environmental_variable} occuring in pathname will be
369: replaced with appropriate values.
371: .seealso: EPSRegisterAll(), EPSRegisterDestroy()
373: M*/
377: PetscErrorCode EPSRegister(char *sname,char *path,char *name,int (*function)(EPS))
378: {
379: PetscErrorCode ierr;
380: char fullname[256];
383: PetscFListConcat(path,name,fullname);
384: PetscFListAdd(&EPSList,sname,fullname,(void (*)(void))function);
385: return(0);
386: }
390: /*@
391: EPSDestroy - Destroys the EPS context.
393: Collective on EPS
395: Input Parameter:
396: . eps - eigensolver context obtained from EPSCreate()
398: Level: beginner
400: .seealso: EPSCreate(), EPSSetUp(), EPSSolve()
401: @*/
402: PetscErrorCode EPSDestroy(EPS eps)
403: {
404: PetscErrorCode ierr;
408: if (--eps->refct > 0) return(0);
410: /* if memory was published with AMS then destroy it */
411: PetscObjectDepublish(eps);
413: STDestroy(eps->OP);
415: if (eps->ops->destroy) {
416: (*eps->ops->destroy)(eps);
417: }
418:
419: if (eps->T) {
420: PetscFree(eps->T);
421: }
422:
423: if (eps->perm) {
424: PetscFree(eps->perm);
425: }
427: if (eps->vec_initial) {
428: VecDestroy(eps->vec_initial);
429: }
431: if (eps->nds > 0) {
432: VecDestroyVecs(eps->DS, eps->nds);
433: }
434:
435: if (eps->DSV) {
436: PetscFree(eps->DSV);
437: }
439: PetscLogObjectDestroy(eps);
440: PetscHeaderDestroy(eps);
441: return(0);
442: }
446: /*@
447: EPSSetST - Associates a spectral transformation object to the
448: eigensolver.
450: Collective on EPS
452: Input Parameters:
453: + eps - eigensolver context obtained from EPSCreate()
454: - st - the spectral transformation object
456: Note:
457: Use EPSGetST() to retrieve the spectral transformation context (for example,
458: to free it at the end of the computations).
460: Level: advanced
462: .seealso: EPSGetST()
463: @*/
464: PetscErrorCode EPSSetST(EPS eps,ST st)
465: {
466: PetscErrorCode ierr;
472: STDestroy(eps->OP);
473: eps->OP = st;
474: PetscObjectReference((PetscObject)eps->OP);
475: return(0);
476: }
480: /*@
481: EPSGetST - Obtain the spectral transformation (ST) object associated
482: to the eigensolver object.
484: Not Collective
486: Input Parameters:
487: . eps - eigensolver context obtained from EPSCreate()
489: Output Parameter:
490: . st - spectral transformation context
492: Level: beginner
494: .seealso: EPSSetST()
495: @*/
496: PetscErrorCode EPSGetST(EPS eps, ST *st)
497: {
500: *st = eps->OP;
501: return(0);
502: }
506: /*@
507: EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
508: eigenvalue problem.
510: Not collective
512: Input Parameter:
513: . eps - the eigenproblem solver context
515: Output Parameter:
516: . is - the answer
518: Level: intermediate
520: @*/
521: PetscErrorCode EPSIsGeneralized(EPS eps,PetscTruth* is)
522: {
523: PetscErrorCode ierr;
524: Mat B;
528: STGetOperators(eps->OP,PETSC_NULL,&B);
529: if( B ) *is = PETSC_TRUE;
530: else *is = PETSC_FALSE;
531: if( eps->setupcalled ) {
532: if( eps->isgeneralized != *is ) {
533: SETERRQ(0,"Warning: Inconsistent EPS state");
534: }
535: }
536: return(0);
537: }
541: /*@
542: EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
543: eigenvalue problem.
545: Not collective
547: Input Parameter:
548: . eps - the eigenproblem solver context
550: Output Parameter:
551: . is - the answer
553: Level: intermediate
555: @*/
556: PetscErrorCode EPSIsHermitian(EPS eps,PetscTruth* is)
557: {
560: if( eps->ishermitian ) *is = PETSC_TRUE;
561: else *is = PETSC_FALSE;
562: return(0);
563: }