Actual source code: itcreate.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: int EPSView(EPS eps,PetscViewer viewer)
39: {
40: char *type, *which;
41: int ierr;
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 basis 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: if (eps->dropvectors) { PetscViewerASCIIPrintf(viewer," computing only eigenvalues\n");}
93: PetscViewerASCIIPrintf(viewer," orthogonalization method: ");
94: switch (eps->orth_type) {
95: case EPS_MGS_ORTH:
96: PetscViewerASCIIPrintf(viewer,"modified Gram-Schmidt\n");
97: break;
98: case EPS_CGS_ORTH:
99: PetscViewerASCIIPrintf(viewer,"classical Gram-Schmidt\n");
100: break;
101: default: SETERRQ(1,"Wrong value of eps->orth_type");
102: }
103: PetscViewerASCIIPrintf(viewer," orthogonalization refinement: ");
104: if (eps->orth_eta == 0.0) {
105: PetscViewerASCIIPrintf(viewer,"never\n");
106: } else if (eps->orth_eta < 0.0) {
107: PetscViewerASCIIPrintf(viewer,"always\n");
108: } else {
109: PetscViewerASCIIPrintf(viewer,"if needed (eta: %f)\n",eps->orth_eta);
110: }
111: PetscViewerASCIIPushTab(viewer);
112: STView(eps->OP,viewer);
113: PetscViewerASCIIPopTab(viewer);
114: } else {
115: if (eps->ops->view) {
116: (*eps->ops->view)(eps,viewer);
117: }
118: STView(eps->OP,viewer);
119: }
120: return(0);
121: }
125: /*@C
126: EPSSetDropEigenvectors - Sets the EPS solver not to compute the
127: eigenvectors. In some methods, this can reduce the number of operations
128: necessary for obtaining the eigenvalues.
130: Collective on KSP
132: Input Parameter:
133: . eps - the eigensolver context
135: Options Database Keys:
136: . -eps_drop_eigenvectors - do not compute eigenvectors
138: Level: advanced
140: .seealso: EPSSetUp(), EPSSolve(), EPSDestroy()
141: @*/
142: int EPSSetDropEigenvectors(EPS eps)
143: {
146: eps->dropvectors = PETSC_TRUE;
147: return(0);
148: }
152: static int EPSPublish_Petsc(PetscObject object)
153: {
154: #if defined(PETSC_HAVE_AMS)
155: EPS v = (EPS) object;
156: int ierr;
157: #endif
158:
161: #if defined(PETSC_HAVE_AMS)
162: /* if it is already published then return */
163: if (v->amem >=0 ) return(0);
165: PetscObjectPublishBaseBegin(object);
166: AMS_Memory_add_field((AMS_Memory)v->amem,"Iteration",&v->its,1,AMS_INT,AMS_READ,
167: AMS_COMMON,AMS_REDUCT_UNDEF);
168: PetscObjectPublishBaseEnd(object);
169: #endif
171: return(0);
172: }
176: /*@C
177: EPSCreate - Creates the default EPS context.
179: Collective on MPI_Comm
181: Input Parameter:
182: . comm - MPI communicator
184: Output Parameter:
185: . eps - location to put the EPS context
187: Note:
188: The default EPS type is EPSPOWER
190: Level: beginner
192: .seealso: EPSSetUp(), EPSSolve(), EPSDestroy(), EPS
193: @*/
194: int EPSCreate(MPI_Comm comm,EPS *outeps)
195: {
196: EPS eps;
197: int ierr;
200: *outeps = 0;
201: PetscHeaderCreate(eps,_p_EPS,struct _EPSOps,EPS_COOKIE,-1,"EPS",comm,EPSDestroy,EPSView);
202: PetscLogObjectCreate(eps);
203: *outeps = eps;
205: eps->bops->publish = EPSPublish_Petsc;
206: eps->ops->setfromoptions = 0;
207: eps->ops->solve = 0;
208: eps->ops->setup = 0;
209: eps->ops->destroy = 0;
210: eps->ops->backtransform = 0;
212: eps->type = -1;
213: eps->max_it = 0;
214: eps->nev = 1;
215: eps->ncv = 0;
216: eps->tol = 0.0;
217: eps->which = EPS_LARGEST_MAGNITUDE;
218: eps->dropvectors = PETSC_FALSE;
219: eps->problem_type = (EPSProblemType)0;
221: eps->vec_initial = 0;
222: eps->V = 0;
223: eps->eigr = 0;
224: eps->eigi = 0;
225: eps->errest = 0;
226: eps->OP = 0;
227: eps->data = 0;
228: eps->nconv = 0;
229: eps->its = 0;
230: eps->perm = PETSC_NULL;
232: eps->nwork = 0;
233: eps->work = 0;
234: eps->isgeneralized = PETSC_FALSE;
235: eps->ishermitian = PETSC_FALSE;
236: eps->setupcalled = 0;
237: eps->reason = EPS_CONVERGED_ITERATING;
239: eps->numbermonitors = 0;
240: eps->numbervmonitors = 0;
242: eps->orthog = EPSClassicalGramSchmidtOrthogonalization;
243: eps->orth_type = EPS_CGS_ORTH;
244: eps->orth_eta = 0.7071;
246: STCreate(comm,&eps->OP);
247: PetscLogObjectParent(eps,eps->OP);
248: PetscPublishAll(eps);
249: return(0);
250: }
251:
254: /*@C
255: EPSSetType - Selects the particular solver to be used in the EPS object.
257: Collective on EPS
259: Input Parameters:
260: + eps - the eigensolver context
261: - type - a known method
263: Options Database Key:
264: . -eps_type <method> - Sets the method; use -help for a list
265: of available methods
266:
267: Notes:
268: See "slepc/include/slepceps.h" for available methods.
270: Normally, it is best to use the EPSSetFromOptions() command and
271: then set the EPS type from the options database rather than by using
272: this routine. Using the options database provides the user with
273: maximum flexibility in evaluating the different available methods.
274: The EPSSetType() routine is provided for those situations where it
275: is necessary to set the iterative solver independently of the command
276: line or options database.
278: Level: intermediate
280: .seealso: STSetType(), EPSType
281: @*/
282: int EPSSetType(EPS eps,EPSType type)
283: {
284: int ierr,(*r)(EPS);
285: PetscTruth match;
291: PetscTypeCompare((PetscObject)eps,type,&match);
292: if (match) return(0);
294: if (eps->data) {
295: /* destroy the old private EPS context */
296: (*eps->ops->destroy)(eps);
297: eps->data = 0;
298: }
299: /* Get the function pointers for the iterative method requested */
300: if (!EPSRegisterAllCalled) {EPSRegisterAll(PETSC_NULL); }
302: PetscFListFind(eps->comm,EPSList,type,(void (**)(void)) &r);
304: if (!r) SETERRQ1(1,"Unknown EPS type given: %s",type);
306: eps->setupcalled = 0;
307: (*r)(eps);
309: PetscObjectChangeTypeName((PetscObject)eps,type);
310: return(0);
311: }
315: /*@C
316: EPSRegisterDestroy - Frees the list of EPS methods that were
317: registered by EPSRegisterDynamic().
319: Not Collective
321: Level: advanced
323: .seealso: EPSRegisterDynamic(), EPSRegisterAll()
324: @*/
325: int EPSRegisterDestroy(void)
326: {
330: if (EPSList) {
331: PetscFListDestroy(&EPSList);
332: EPSList = 0;
333: }
334: EPSRegisterAllCalled = PETSC_FALSE;
335: return(0);
336: }
340: /*@C
341: EPSGetType - Gets the EPS type as a string from the EPS object.
343: Not Collective
345: Input Parameter:
346: . eps - the eigensolver context
348: Output Parameter:
349: . name - name of EPS method
351: Level: intermediate
353: .seealso: EPSSetType()
354: @*/
355: int EPSGetType(EPS eps,EPSType *type)
356: {
359: *type = eps->type_name;
360: return(0);
361: }
365: /*@
366: EPSSetFromOptions - Sets EPS options from the options database.
367: This routine must be called before EPSSetUp() if the user is to be
368: allowed to set the solver type.
370: Collective on EPS
372: Input Parameters:
373: . eps - the eigensolver context
375: Notes:
376: To see all options, run your program with the -help option.
378: Level: beginner
380: .seealso:
381: @*/
382: int EPSSetFromOptions(EPS eps)
383: {
384: int ierr;
385: char type[256];
386: PetscTruth flg;
387: const char *orth_list[2] = { "mgs" , "cgs" };
388: const char *ref_list[3] = { "never" , "ifneeded", "always" };
389: EPSOrthogonalizationType orth_type;
390: EPSOrthogonalizationRefinementType ref_type;
391: PetscReal eta;
395: if (!EPSRegisterAllCalled) {EPSRegisterAll(PETSC_NULL);}
396: PetscOptionsBegin(eps->comm,eps->prefix,"Eigenproblem Solver (EPS) Options","EPS");
397: PetscOptionsList("-eps_type","Eigenproblem Solver method","EPSSetType",EPSList,(char*)(eps->type_name?eps->type_name:EPSPOWER),type,256,&flg);
398: if (flg) {
399: EPSSetType(eps,type);
400: }
401: /*
402: Set the type if it was never set.
403: */
404: if (!eps->type_name) {
405: EPSSetType(eps,EPSPOWER);
406: }
408: PetscOptionsLogicalGroupBegin("-eps_hermitian","hermitian eigenvalue problem","EPSSetProblemType",&flg);
409: if (flg) {EPSSetProblemType(eps,EPS_HEP);}
410: PetscOptionsLogicalGroup("-eps_gen_hermitian","generalized hermitian eigenvalue problem","EPSSetProblemType",&flg);
411: if (flg) {EPSSetProblemType(eps,EPS_GHEP);}
412: PetscOptionsLogicalGroup("-eps_non_hermitian","non-hermitian eigenvalue problem","EPSSetProblemType",&flg);
413: if (flg) {EPSSetProblemType(eps,EPS_NHEP);}
414: PetscOptionsLogicalGroupEnd("-eps_gen_non_hermitian","generalized non-hermitian eigenvalue problem","EPSSetProblemType",&flg);
415: if (flg) {EPSSetProblemType(eps,EPS_GNHEP);}
417: orth_type = eps->orth_type;
418: PetscOptionsEList("-eps_orthog_type","Orthogonalization method","EPSSetOrthogonalization",orth_list,2,orth_list[orth_type],(int*)&orth_type,&flg);
419: ref_type = EPS_ORTH_REFINE_IFNEEDED;
420: PetscOptionsEList("-eps_orthog_refinement","Iterative refinement mode during orthogonalization","EPSSetOrthogonalizationRefinement",ref_list,3,ref_list[ref_type],(int*)&ref_type,&flg);
421: eta = eps->orth_eta;
422: PetscOptionsReal("-eps_orthog_eta","Parameter of iterative refinement during orthogonalization","EPSSetOrthogonalizationRefinement",eta,&eta,PETSC_NULL);
423: EPSSetOrthogonalization(eps,orth_type,ref_type,eta);
425: PetscOptionsInt("-eps_max_it","Maximum number of iterations","EPSSetTolerances",eps->max_it,&eps->max_it,PETSC_NULL);
426: PetscOptionsReal("-eps_tol","Tolerance","KSPSetTolerances",eps->tol,&eps->tol,PETSC_NULL);
427: PetscOptionsInt("-eps_nev","Number of eigenvalues to compute","EPSSetDimensions",eps->nev,&eps->nev,&flg);
428: if( eps->nev<1 ) SETERRQ(1,"Illegal value for option -eps_nev. Must be > 0");
429: PetscOptionsInt("-eps_ncv","Number of basis vectors","EPSSetDimensions",eps->ncv,&eps->ncv,&flg);
430: if( flg && eps->ncv<1 ) SETERRQ(1,"Illegal value for option -eps_ncv. Must be > 0");
432: PetscOptionsName("-eps_drop_eigenvectors","Do not compute eigenvectors","EPSSetDropEigenvectors",&flg);
433: if (flg) {
434: EPSSetDropEigenvectors(eps);
435: }
437: /* -----------------------------------------------------------------------*/
438: /*
439: Cancels all monitors hardwired into code before call to EPSSetFromOptions()
440: */
441: PetscOptionsName("-eps_cancelmonitors","Remove any hardwired monitor routines","EPSClearMonitor",&flg);
442: if (flg) {
443: EPSClearMonitor(eps);
444: }
445: /*
446: Prints error estimates at each iteration
447: */
448: PetscOptionsName("-eps_monitor","Monitor error estimates","EPSSetMonitor",&flg);
449: if (flg) {
450: EPSSetMonitor(eps,EPSDefaultEstimatesMonitor,PETSC_NULL);
451: }
452: /*
453: Prints approximate eigenvalues at each iteration
454: */
455: PetscOptionsName("-eps_monitor_values","Monitor approximate eigenvalues","EPSSetValuesMonitor",&flg);
456: if (flg) {
457: EPSSetValuesMonitor(eps,EPSDefaultValuesMonitor,PETSC_NULL);
458: }
459: /* -----------------------------------------------------------------------*/
460: PetscOptionsLogicalGroupBegin("-eps_largest_magnitude","compute largest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg);
461: if (flg) {EPSSetWhichEigenpairs(eps,EPS_LARGEST_MAGNITUDE);}
462: PetscOptionsLogicalGroup("-eps_smallest_magnitude","compute smallest eigenvalues in magnitude","EPSSetWhichEigenpairs",&flg);
463: if (flg) {EPSSetWhichEigenpairs(eps,EPS_SMALLEST_MAGNITUDE);}
464: PetscOptionsLogicalGroup("-eps_largest_real","compute largest real parts","EPSSetWhichEigenpairs",&flg);
465: if (flg) {EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL);}
466: PetscOptionsLogicalGroup("-eps_smallest_real","compute smallest real parts","EPSSetWhichEigenpairs",&flg);
467: if (flg) {EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL);}
468: PetscOptionsLogicalGroup("-eps_largest_imaginary","compute largest imaginary parts","EPSSetWhichEigenpairs",&flg);
469: if (flg) {EPSSetWhichEigenpairs(eps,EPS_LARGEST_IMAGINARY);}
470: PetscOptionsLogicalGroupEnd("-eps_smallest_imaginary","compute smallest imaginary parts","EPSSetWhichEigenpairs",&flg);
471: if (flg) {EPSSetWhichEigenpairs(eps,EPS_SMALLEST_IMAGINARY);}
473: PetscOptionsName("-eps_view","Print detailed information on solver used","EPSView",0);
474: PetscOptionsName("-eps_view_binary","Save the matrices associated to the eigenproblem","EPSSetFromOptions",0);
475: PetscOptionsName("-eps_plot_eigs","Make a plot of the computed eigenvalues","EPSSolve",0);
477: if (eps->ops->setfromoptions) {
478: (*eps->ops->setfromoptions)(eps);
479: }
480: PetscOptionsEnd();
482: STSetFromOptions(eps->OP);
483: return(0);
484: }
486: /*MC
487: EPSRegisterDynamic - Adds a method to the eigenproblem solver package.
489: Synopsis:
490: EPSRegisterDynamic(char *name_solver,char *path,char *name_create,int (*routine_create)(EPS))
492: Not Collective
494: Input Parameters:
495: + name_solver - name of a new user-defined solver
496: . path - path (either absolute or relative) the library containing this solver
497: . name_create - name of routine to create the solver context
498: - routine_create - routine to create the solver context
500: Notes:
501: EPSRegisterDynamic() may be called multiple times to add several user-defined solvers.
503: If dynamic libraries are used, then the fourth input argument (routine_create)
504: is ignored.
506: Sample usage:
507: .vb
508: EPSRegisterDynamic("my_solver",/home/username/my_lib/lib/libO/solaris/mylib.a,
509: "MySolverCreate",MySolverCreate);
510: .ve
512: Then, your solver can be chosen with the procedural interface via
513: $ EPSSetType(eps,"my_solver")
514: or at runtime via the option
515: $ -eps_type my_solver
517: Level: advanced
519: Environmental variables such as ${PETSC_ARCH}, ${SLEPC_DIR}, ${BOPT},
520: and others of the form ${any_environmental_variable} occuring in pathname will be
521: replaced with appropriate values.
523: .seealso: EPSRegisterAll(), EPSRegisterDestroy()
525: M*/
529: int EPSRegister(char *sname,char *path,char *name,int (*function)(EPS))
530: {
531: int ierr;
532: char fullname[256];
535: PetscFListConcat(path,name,fullname);
536: PetscFListAdd(&EPSList,sname,fullname,(void (*)(void))function);
537: return(0);
538: }