Actual source code: itfunc.c
1: /*
2: Interface EPS routines that the user calls.
3: */
5: #include src/eps/epsimpl.h
6: #include slepcblaslapack.h
10: /*@
11: EPSSetDefaults - Sets the default values of the internal parameters
12: of the eigensolver. Some of these parameters are set by a method-specific
13: routine.
15: Collective on EPS
17: Input Parameter:
18: . eps - eigenproblem solver context
20: Level: developer
22: .seealso: EPSSetType(), EPSSetUp()
23: @*/
24: int EPSSetDefaults(EPS eps)
25: {
26: int ierr;
27: Mat A,B;
28: PetscTruth Ah,Bh;
29:
33: /* Set default solver type */
34: if (!eps->type_name) {
35: EPSSetType(eps,EPSPOWER);
36: }
38: /* Set default problem type */
39: if (!eps->problem_type) {
40: STGetOperators(eps->OP,&A,&B);
41: SlepcIsHermitian(A,&Ah);
42: if (B==PETSC_NULL) {
43: if (Ah) { EPSSetProblemType(eps,EPS_HEP); }
44: else { EPSSetProblemType(eps,EPS_NHEP); }
45: }
46: else {
47: SlepcIsHermitian(B,&Bh);
48: if (Ah && Bh) { EPSSetProblemType(eps,EPS_GHEP); }
49: else { EPSSetProblemType(eps,EPS_GNHEP); }
50: }
51: }
53: if (eps->ops->setdefaults) {
54: (*eps->ops->setdefaults)(eps);
55: }
56: return(0);
57: }
61: /*@
62: EPSSetUp - Sets up all the internal data structures necessary for the
63: execution of the eigensolver. Then calls STSetUp() for any set-up
64: operations associated to the ST object.
66: Collective on EPS
68: Input Parameter:
69: . eps - eigenproblem solver context
71: Level: advanced
73: Notes:
74: This function need not be called explicitly in most cases, since EPSSolve()
75: calls it. It can be useful when one wants to measure the set-up time
76: separately from the solve time.
78: This function sets a random initial vector if none has been provided.
80: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp()
81: @*/
82: int EPSSetUp(EPS eps)
83: {
84: int n, m, i, ierr, nloc, nev, ncv;
85: PetscScalar *pV;
86: Vec v0;
87: Mat A;
88:
92: /* reset the convergence flag from the previous solves */
93: eps->reason = EPS_CONVERGED_ITERATING;
95: /* Check if the EPS initial vector has been set */
96: EPSGetInitialVector(eps,&v0);
97: if (!v0) {
98: STGetOperators(eps->OP,&A,PETSC_NULL);
99: MatGetLocalSize(A,&n,&m);
100: VecCreate(eps->comm,&v0);
101: VecSetSizes(v0,m,PETSC_DECIDE);
102: VecSetFromOptions(v0);
103: SlepcVecSetRandom(v0);
104: eps->vec_initial = v0;
105: }
106: STSetVector(eps->OP,v0);
108: EPSSetDefaults(eps);
109:
110: nev = eps->nev;
111: ncv = eps->ncv;
112: if (!eps->eigr){
113: PetscMalloc(ncv*sizeof(PetscScalar),&eps->eigr);
114: }
115: if (!eps->eigi){
116: PetscMalloc(ncv*sizeof(PetscScalar),&eps->eigi);
117: }
118: if (!eps->errest){
119: PetscMalloc(ncv*sizeof(PetscReal),&eps->errest);
120: }
121: if (!eps->V){
122: VecGetLocalSize(eps->vec_initial,&nloc);
123: PetscMalloc(ncv*sizeof(Vec),&eps->V);
124: PetscMalloc(ncv*nloc*sizeof(PetscScalar),&pV);
125: for (i=0;i<ncv;i++) {
126: VecCreateMPIWithArray(eps->comm,nloc,PETSC_DECIDE,pV+i*nloc,&eps->V[i]);
127: }
128: }
130: if (eps->setupcalled) return(0);
131: PetscLogEventBegin(EPS_SetUp,eps,eps->V[0],eps->V[0],0);
132: eps->setupcalled = 1;
133: (*eps->ops->setup)(eps);
134: STSetUp(eps->OP);
135: PetscLogEventEnd(EPS_SetUp,eps,eps->V[0],eps->V[0],0);
136: return(0);
137: }
141: /*@
142: EPSSolve - Solves the eigensystem.
144: Collective on EPS
146: Input Parameter:
147: . eps - eigensolver context obtained from EPSCreate()
149: Options Database:
150: + -eps_view - print information about the solver used
151: . -eps_view_binary - save the matrices to the default binary file
152: - -eps_plot_eigs - plot computed eigenvalues
154: Level: beginner
156: .seealso: EPSCreate(), EPSSetUp(), EPSDestroy(), EPSSetTolerances()
157: @*/
158: int EPSSolve(EPS eps)
159: {
160: int i,ierr;
161: PetscReal re,im;
162: PetscTruth flg;
163: PetscViewer viewer;
164: PetscDraw draw;
165: PetscDrawSP drawsp;
170: PetscOptionsHasName(eps->prefix,"-eps_view_binary",&flg);
171: if (flg) {
172: Mat A,B;
173: STGetOperators(eps->OP,&A,&B);
174: MatView(A,PETSC_VIEWER_BINARY_(eps->comm));
175: if (B) MatView(B,PETSC_VIEWER_BINARY_(eps->comm));
176: }
178: if (!eps->setupcalled){ EPSSetUp(eps); }
179: STResetNumberLinearIterations(eps->OP);
180: PetscLogEventBegin(EPS_Solve,eps,eps->V[0],eps->V[0],0);
181: STPreSolve(eps->OP,eps);
182: (*eps->ops->solve)(eps);
183: STPostSolve(eps->OP,eps);
184: PetscLogEventEnd(EPS_Solve,eps,eps->V[0],eps->V[0],0);
185: if (!eps->reason) {
186: SETERRQ(1,"Internal error, solver returned without setting converged reason");
187: }
189: /* Map eigenvalues back to the original problem, necessary in some
190: * spectral transformations */
191: (*eps->ops->backtransform)(eps);
193: #ifndef PETSC_USE_COMPLEX
194: /* reorder conjugate eigenvalues (positive imaginary first) */
195: for (i=0; i<eps->nconv-1; i++) {
196: PetscScalar minus = -1.0;
197: if (eps->eigi[i] != 0) {
198: if (eps->eigi[i] < 0) {
199: eps->eigi[i] = -eps->eigi[i];
200: eps->eigi[i+1] = -eps->eigi[i+1];
201: VecScale(&minus, eps->V[i+1]);
202: }
203: i++;
204: }
205: }
206: #endif
208: /* sort eigenvalues according to eps->which parameter */
209: if (eps->nconv > 0) {
210: PetscMalloc(sizeof(int)*eps->nconv, &eps->perm);
211: EPSSortEigenvalues(eps->nconv, eps->eigr, eps->eigi, eps->which, eps->nconv, eps->perm);
212: }
214: PetscOptionsHasName(eps->prefix,"-eps_view",&flg);
215: if (flg && !PetscPreLoadingOn) { EPSView(eps,PETSC_VIEWER_STDOUT_WORLD); }
217: PetscOptionsHasName(eps->prefix,"-eps_plot_eigs",&flg);
218: if (flg) {
219: PetscViewerDrawOpen(PETSC_COMM_SELF,0,"Computed Eigenvalues",
220: PETSC_DECIDE,PETSC_DECIDE,300,300,&viewer);
221: PetscViewerDrawGetDraw(viewer,0,&draw);
222: PetscDrawSPCreate(draw,1,&drawsp);
223: for( i=0; i<eps->nconv; i++ ) {
224: #if defined(PETSC_USE_COMPLEX)
225: re = PetscRealPart(eps->eigr[i]);
226: im = PetscImaginaryPart(eps->eigi[i]);
227: #else
228: re = eps->eigr[i];
229: im = eps->eigi[i];
230: #endif
231: PetscDrawSPAddPoint(drawsp,&re,&im);
232: }
233: PetscDrawSPDraw(drawsp);
234: PetscDrawSPDestroy(drawsp);
235: PetscViewerDestroy(viewer);
236: }
238: return(0);
239: }
243: /*@
244: EPSDestroy - Destroys the EPS context.
246: Collective on EPS
248: Input Parameter:
249: . eps - eigensolver context obtained from EPSCreate()
251: Level: beginner
253: .seealso: EPSCreate(), EPSSetUp(), EPSSolve()
254: @*/
255: int EPSDestroy(EPS eps)
256: {
257: int i,ierr;
258: PetscScalar *pV;
262: if (--eps->refct > 0) return(0);
264: /* if memory was published with AMS then destroy it */
265: PetscObjectDepublish(eps);
267: STDestroy(eps->OP);
268: if (eps->ops->destroy) {
269: (*eps->ops->destroy)(eps);
270: }
272: if (eps->eigr){ PetscFree(eps->eigr); }
273: if (eps->eigi){ PetscFree(eps->eigi); }
274: if (eps->errest){ PetscFree(eps->errest); }
275: if (eps->V){
276: VecGetArray(eps->V[0],&pV);
277: for (i=0;i<eps->ncv;i++) {
278: VecDestroy(eps->V[i]);
279: }
280: PetscFree(pV);
281: PetscFree(eps->V);
282: }
283: if (eps->vec_initial) {
284: VecDestroy(eps->vec_initial);
285: }
286: if (eps->perm) {
287: PetscFree(eps->perm);
288: }
290: PetscLogObjectDestroy(eps);
291: PetscHeaderDestroy(eps);
292: return(0);
293: }
297: /*@
298: EPSGetTolerances - Gets the tolerance and maximum
299: iteration count used by the default EPS convergence tests.
301: Not Collective
303: Input Parameter:
304: . eps - the eigensolver context
305:
306: Output Parameters:
307: + tol - the convergence tolerance
308: - maxits - maximum number of iterations
310: Notes:
311: The user can specify PETSC_NULL for any parameter that is not needed.
313: Level: intermediate
315: .seealso: EPSSetTolerances()
316: @*/
317: int EPSGetTolerances(EPS eps,PetscReal *tol,int *maxits)
318: {
321: if (tol) *tol = eps->tol;
322: if (maxits) *maxits = eps->max_it;
323: return(0);
324: }
328: /*@
329: EPSSetTolerances - Sets the tolerance and maximum
330: iteration count used by the default EPS convergence testers.
332: Collective on EPS
334: Input Parameters:
335: + eps - the eigensolver context
336: . tol - the convergence tolerance
337: - maxits - maximum number of iterations to use
339: Options Database Keys:
340: + -eps_tol <tol> - Sets the convergence tolerance
341: - -eps_max_it <maxits> - Sets the maximum number of iterations allowed
343: Notes:
344: Use PETSC_DEFAULT to retain the default value of any of the tolerances.
346: Level: intermediate
348: .seealso: EPSGetTolerances()
349: @*/
350: int EPSSetTolerances(EPS eps,PetscReal tol,int maxits)
351: {
354: if (tol != PETSC_DEFAULT) eps->tol = tol;
355: if (maxits != PETSC_DEFAULT) eps->max_it = maxits;
356: return(0);
357: }
361: /*@
362: EPSGetDimensions - Gets the number of eigenvalues to compute
363: and the dimension of the subspace.
365: Not Collective
367: Input Parameter:
368: . eps - the eigensolver context
369:
370: Output Parameters:
371: + nev - number of eigenvalues to compute
372: - ncv - the maximum dimension of the subspace to be used by the solver
374: Notes:
375: The user can specify PETSC_NULL for any parameter that is not needed.
377: Level: intermediate
379: .seealso: EPSSetDimensions()
380: @*/
381: int EPSGetDimensions(EPS eps,int *nev,int *ncv)
382: {
385: if( nev ) *nev = eps->nev;
386: if( ncv ) *ncv = eps->ncv;
387: return(0);
388: }
392: /*@
393: EPSSetDimensions - Sets the number of eigenvalues to compute
394: and the dimension of the subspace.
396: Collective on EPS
398: Input Parameters:
399: + eps - the eigensolver context
400: . nev - number of eigenvalues to compute
401: - ncv - the maximum dimension of the subspace to be used by the solver
403: Options Database Keys:
404: + -eps_nev <nev> - Sets the number of eigenvalues
405: - -eps_ncv <ncv> - Sets the dimension of the subspace
407: Notes:
408: Use PETSC_DEFAULT to retain the previous value of any parameter.
410: Use PETSC_DECIDE for ncv to assign a reasonably good value, which is
411: dependent on the solution method.
413: Level: intermediate
415: .seealso: EPSGetDimensions()
416: @*/
417: int EPSSetDimensions(EPS eps,int nev,int ncv)
418: {
422: if( nev != PETSC_DEFAULT ) {
423: if (nev<1) SETERRQ(1,"Illegal value of nev. Must be > 0");
424: eps->nev = nev;
425: }
426: if( ncv == PETSC_DECIDE ) eps->ncv = 0;
427: else if( ncv != PETSC_DEFAULT ) {
428: if (ncv<1) SETERRQ(1,"Illegal value of ncv. Must be > 0");
429: eps->ncv = ncv;
430: }
431: return(0);
432: }
436: /*@
437: EPSGetConverged - Gets the number of converged eigenpairs.
439: Not Collective
441: Input Parameter:
442: . eps - the eigensolver context
443:
444: Output Parameter:
445: . nconv - number of converged eigenpairs
447: Note:
448: This function should be called after EPSSolve() has finished.
450: Level: beginner
452: .seealso: EPSSetDimensions()
453: @*/
454: int EPSGetConverged(EPS eps,int *nconv)
455: {
458: if (nconv) *nconv = eps->nconv;
459: return(0);
460: }
464: /*@
465: EPSGetEigenpair - Gets the i-th solution of the eigenproblem
466: as computed by EPSSolve(). The solution consists in both the eigenvalue
467: and the eigenvector (if available).
469: Not Collective
471: Input Parameters:
472: + eps - eigensolver context
473: - i - index of the solution
475: Output Parameters:
476: + eigr - real part of eigenvalue
477: . eigi - imaginary part of eigenvalue
478: . Vr - real part of eigenvector
479: - Vi - imaginary part of eigenvector
481: Notes:
482: If the eigenvalue is real, then eigi and Vi are set to zero. In the
483: complex case (e.g. with BOPT=O_complex) the eigenvalue is stored
484: directly in eigr (eigi is set to zero) and the eigenvector Vr (Vi is
485: set to zero).
487: The index i should be a value between 0 and nconv (see EPSGetConverged()).
488: Eigenpairs are indexed according to the ordering criterion established
489: with EPSSetWhichEigenpairs().
491: Level: beginner
493: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
494: @*/
495: int EPSGetEigenpair(EPS eps, int i, PetscScalar *eigr, PetscScalar *eigi, Vec Vr, Vec Vi)
496: {
497: int ierr, k;
498: PetscScalar zero = 0.0, minus = -1.0;
502: if (!eps->eigr || !eps->eigi || !eps->V) {
503: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
504: }
505: if (i<0 || i>=eps->nconv) {
506: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
507: }
508: if (eps->dropvectors && (Vr || Vi) ) {
509: SETERRQ(1, "Eigenvectors are not available");
510: }
512: if (!eps->perm) k = i;
513: else k = eps->perm[i];
514: #ifdef PETSC_USE_COMPLEX
515: if (eigr) *eigr = eps->eigr[k];
516: if (eigi) *eigi = 0;
517: if (Vr) { VecCopy(eps->V[k], Vr); }
518: if (Vi) { VecSet(&zero, Vi); }
519: #else
520: if (eigr) *eigr = eps->eigr[k];
521: if (eigi) *eigi = eps->eigi[k];
522: if (eps->eigi[k] > 0) { /* first value of conjugate pair */
523: if (Vr) { VecCopy(eps->V[k], Vr); }
524: if (Vi) { VecCopy(eps->V[k+1], Vi); }
525: } else if (eps->eigi[k] < 0) { /* second value of conjugate pair */
526: if (Vr) { VecCopy(eps->V[k-1], Vr); }
527: if (Vi) {
528: VecCopy(eps->V[k], Vi);
529: VecScale(&minus, Vi);
530: }
531: } else { /* real eigenvalue */
532: if (Vr) { VecCopy(eps->V[k], Vr); }
533: if (Vi) { VecSet(&zero, Vi); }
534: }
535: #endif
536:
537: return(0);
538: }
542: /*@
543: EPSGetErrorEstimate - Returns the error bound associated to the i-th
544: approximate eigenpair.
546: Not Collective
548: Input Parameter:
549: + eps - eigensolver context
550: - i - index of eigenpair
552: Output Parameter:
553: . errest - the error estimate
555: Level: advanced
557: .seealso: EPSComputeRelativeError()
558: @*/
559: int EPSGetErrorEstimate(EPS eps, int i, PetscReal *errest)
560: {
563: if (!eps->eigr || !eps->eigi) {
564: SETERRQ(PETSC_ERR_ARG_WRONGSTATE, "EPSSolve must be called first");
565: }
566: if (i<0 || i>=eps->nconv) {
567: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE, "Argument 2 out of range");
568: }
569: if (eps->perm) i = eps->perm[i];
570: if (errest) *errest = eps->errest[i];
571: return(0);
572: }
576: /*@
577: EPSSetST - Associates a spectral transformation object to the
578: eigensolver.
580: Collective on EPS
582: Input Parameters:
583: + eps - eigensolver context obtained from EPSCreate()
584: - st - the spectral transformation object
586: Note:
587: Use EPSGetST() to retrieve the spectral transformation context (for example,
588: to free it at the end of the computations).
590: Level: advanced
592: .seealso: EPSGetST()
593: @*/
594: int EPSSetST(EPS eps,ST st)
595: {
602: STDestroy(eps->OP);
603: eps->OP = st;
604: PetscObjectReference((PetscObject)eps->OP);
605: return(0);
606: }
610: /*@
611: EPSGetST - Obtain the spectral transformation (ST) object associated
612: to the eigensolver object.
614: Not Collective
616: Input Parameters:
617: . eps - eigensolver context obtained from EPSCreate()
619: Output Parameter:
620: . st - spectral transformation context
622: Level: beginner
624: .seealso: EPSSetST()
625: @*/
626: int EPSGetST(EPS eps, ST *st)
627: {
630: *st = eps->OP;
631: return(0);
632: }
636: /*@C
637: EPSSetMonitor - Sets an ADDITIONAL function to be called at every
638: iteration to monitor the error estimates for each requested eigenpair.
639:
640: Collective on EPS
642: Input Parameters:
643: + eps - eigensolver context obtained from EPSCreate()
644: . monitor - pointer to function (if this is PETSC_NULL, it turns off monitoring
645: - mctx - [optional] context for private data for the
646: monitor routine (use PETSC_NULL if no context is desired)
648: Calling Sequence of monitor:
649: $ monitor (EPS eps, int its, int nconv, PetscReal* errest, int nest, void *mctx)
651: + eps - eigensolver context obtained from EPSCreate()
652: . its - iteration number
653: . nconv - number of converged eigenpairs
654: . errest - error estimates for each eigenpair
655: . nest - number of error estimates
656: - mctx - optional monitoring context, as set by EPSSetMonitor()
658: Options Database Keys:
659: + -eps_monitor - print error estimates at each iteration
660: - -eps_cancelmonitors - cancels all monitors that have been hardwired into
661: a code by calls to EPSetMonitor(), but does not cancel those set via
662: the options database.
664: Notes:
665: Several different monitoring routines may be set by calling
666: EPSSetMonitor() multiple times; all will be called in the
667: order in which they were set.
669: Level: intermediate
671: .seealso: EPSDefaultEstimatesMonitor(), EPSClearMonitor()
672: @*/
673: int EPSSetMonitor(EPS eps, int (*monitor)(EPS,int,int,PetscReal*,int,void*), void *mctx)
674: {
677: if (eps->numbermonitors >= MAXEPSMONITORS) {
678: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many EPS monitors set");
679: }
680: eps->monitor[eps->numbermonitors] = monitor;
681: eps->monitorcontext[eps->numbermonitors++] = (void*)mctx;
682: return(0);
683: }
687: /*@C
688: EPSSetValuesMonitor - Sets an ADDITIONAL function to be called at every
689: iteration to monitor the value of approximate eigenvalues.
690:
691: Collective on EPS
693: Input Parameters:
694: + eps - eigensolver context obtained from EPSCreate()
695: . monitor - pointer to function (if this is PETSC_NULL, it turns off monitoring
696: - mctx - [optional] context for private data for the
697: monitor routine (use PETSC_NULL if no context is desired)
699: Calling Sequence of monitor:
700: $ monitor (EPS eps, int its, int nconv, PetscScalar* kr, PetscScalar* ki, int nest, void *mctx)
702: + eps - eigensolver context obtained from EPSCreate()
703: . its - iteration number
704: . nconv - number of converged eigenpairs
705: . kr - real part of each eigenvalue
706: . ki - imaginary part of each eigenvalue
707: . nest - number of error estimates
708: - mctx - optional monitoring context, as set by EPSSetMonitor()
710: Options Database Keys:
711: + -eps_monitor_values - print eigenvalue approximations at each iteration
712: - -eps_cancelmonitors - cancels all monitors that have been hardwired into
713: a code by calls to EPSetValuesMonitor(), but does not cancel those set
714: via the options database.
716: Notes:
717: Several different monitoring routines may be set by calling
718: EPSSetValuesMonitor() multiple times; all will be called in the
719: order in which they were set.
721: Level: intermediate
723: .seealso: EPSDefaultValuesMonitor(), EPSClearMonitor()
724: @*/
725: int EPSSetValuesMonitor(EPS eps, int (*monitor)(EPS,int,int,PetscScalar*,PetscScalar*,int,void*), void *mctx)
726: {
729: if (eps->numbervmonitors >= MAXEPSMONITORS) {
730: SETERRQ(PETSC_ERR_ARG_OUTOFRANGE,"Too many EPS values monitors set");
731: }
732: eps->vmonitor[eps->numbervmonitors] = monitor;
733: eps->vmonitorcontext[eps->numbervmonitors++] = (void*)mctx;
734: return(0);
735: }
739: /*@C
740: EPSClearMonitor - Clears all monitors for an EPS object.
742: Collective on EPS
744: Input Parameters:
745: . eps - eigensolver context obtained from EPSCreate()
747: Options Database Key:
748: . -eps_cancelmonitors - Cancels all monitors that have been hardwired
749: into a code by calls to EPSSetMonitor() or EPSSetValuesMonitor(),
750: but does not cancel those set via the options database.
752: Level: intermediate
754: .seealso: EPSSetMonitor(), EPSSetValuesMonitor()
755: @*/
756: int EPSClearMonitor(EPS eps)
757: {
760: eps->numbermonitors = 0;
761: eps->numbervmonitors = 0;
762: return(0);
763: }
767: /*@C
768: EPSGetMonitorContext - Gets the estimates monitor context, as set by
769: EPSSetMonitor() for the FIRST monitor only.
771: Not Collective
773: Input Parameter:
774: . eps - eigensolver context obtained from EPSCreate()
776: Output Parameter:
777: . ctx - monitor context
779: Level: intermediate
781: .seealso: EPSSetMonitor(), EPSDefaultEstimatesMonitor()
782: @*/
783: int EPSGetMonitorContext(EPS eps, void **ctx)
784: {
787: *ctx = (eps->monitorcontext[0]);
788: return(0);
789: }
793: /*@C
794: EPSGetValuesMonitorContext - Gets the values monitor context, as set by
795: EPSSetValuesMonitor() for the FIRST monitor only.
797: Not Collective
799: Input Parameter:
800: . eps - eigensolver context obtained from EPSCreate()
802: Output Parameter:
803: . ctx - monitor context
805: Level: intermediate
807: .seealso: EPSSetValuesMonitor(), EPSDefaultValuesMonitor()
808: @*/
809: int EPSGetValuesMonitorContext(EPS eps, void **ctx)
810: {
813: *ctx = (eps->vmonitorcontext[0]);
814: return(0);
815: }
819: /*@
820: EPSSetInitialVector - Sets the initial vector from which the
821: eigensolver starts to iterate.
823: Collective on EPS and Vec
825: Input Parameters:
826: + eps - the eigensolver context
827: - vec - the vector
829: Level: intermediate
831: .seealso: EPSGetInitialVector()
833: @*/
834: int EPSSetInitialVector(EPS eps,Vec vec)
835: {
837:
842: if (eps->vec_initial) {
843: VecDestroy(eps->vec_initial);
844: }
845: eps->vec_initial = vec;
846: PetscObjectReference((PetscObject)eps->vec_initial);
847: return(0);
848: }
852: /*@
853: EPSGetInitialVector - Gets the initial vector associated with the
854: eigensolver; if the vector was not set it will return a 0 pointer or
855: a vector randomly generated by EPSSetUp().
857: Not collective, but vector is shared by all processors that share the EPS
859: Input Parameter:
860: . eps - the eigensolver context
862: Output Parameter:
863: . vec - the vector
865: Level: intermediate
867: .seealso: EPSSetInitialVector()
869: @*/
870: int EPSGetInitialVector(EPS eps,Vec *vec)
871: {
874: *vec = eps->vec_initial;
875: return(0);
876: }
880: /*@
881: EPSSetWhichEigenpairs - Specifies which portion of the spectrum is
882: to be sought.
884: Collective on EPS
886: Input Parameter:
887: . eps - eigensolver context obtained from EPSCreate()
889: Output Parameter:
890: . which - the portion of the spectrum to be sought
892: Possible values:
893: The parameter 'which' can have one of these values:
894:
895: + EPS_LARGEST_MAGNITUDE - largest eigenvalues in magnitude (default)
896: . EPS_SMALLEST_MAGNITUDE - smallest eigenvalues in magnitude
897: . EPS_LARGEST_REAL - largest real parts
898: . EPS_SMALLEST_REAL - smallest real parts
899: . EPS_LARGEST_IMAGINARY - largest imaginary parts
900: - EPS_SMALLEST_IMAGINARY - smallest imaginary parts
902: Options Database Keys:
903: + -eps_largest_magnitude - Sets largest eigenvalues in magnitude
904: . -eps_smallest_magnitude - Sets smallest eigenvalues in magnitude
905: . -eps_largest_real - Sets largest real parts
906: . -eps_smallest_real - Sets smallest real parts
907: . -eps_largest_imaginary - Sets largest imaginary parts in magnitude
908: - -eps_smallest_imaginary - Sets smallest imaginary parts in magnitude
910: Notes:
911: Not all eigensolvers implemented in EPS account for all the possible values
912: stated above. Also, some values make sense only for certain types of
913: problems. If SLEPc is compiled for real numbers EPS_LARGEST_IMAGINARY
914: and EPS_SMALLEST_IMAGINARY use the absolute value of the imaginary part
915: for eigenvalue selection.
916:
917: Level: intermediate
919: .seealso: EPSGetWhichEigenpairs(), EPSSortEigenvalues()
920: @*/
921: int EPSSetWhichEigenpairs(EPS eps,EPSWhich which)
922: {
925: eps->which = which;
926: return(0);
927: }
931: /*@
932: EPSGetWhichEigenpairs - Returns which portion of the spectrum is to be
933: sought.
935: Not Collective
937: Input Parameter:
938: . eps - eigensolver context obtained from EPSCreate()
940: Output Parameter:
941: . which - the portion of the spectrum to be sought
943: Notes:
944: See EPSSetWhichEigenpairs() for possible values of which
946: Level: intermediate
948: .seealso: EPSSetWhichEigenpairs()
949: @*/
950: int EPSGetWhichEigenpairs(EPS eps,EPSWhich *which)
951: {
954: *which = eps->which;
955: return(0);
956: }
960: /*@
961: EPSComputeExplicitOperator - Computes the explicit operator associated
962: to the eigenvalue problem with the specified spectral transformation.
964: Collective on EPS
966: Input Parameter:
967: . eps - the eigenvalue solver context
969: Output Parameter:
970: . mat - the explicit operator
972: Notes:
973: This routine builds a matrix containing the explicit operator. For
974: example, in generalized problems with shift-and-invert spectral
975: transformation the result would be matrix (A - s B)^-1 B.
976:
977: This computation is done by applying the operator to columns of the
978: identity matrix.
980: Currently, this routine uses a dense matrix format when 1 processor
981: is used and a sparse format otherwise. This routine is costly in general,
982: and is recommended for use only with relatively small systems.
984: Level: advanced
986: .seealso: STApply()
987: @*/
988: int EPSComputeExplicitOperator(EPS eps,Mat *mat)
989: {
990: Vec in,out;
991: int ierr,i,M,m,size,*rows,start,end;
992: MPI_Comm comm;
993: PetscScalar *array,zero = 0.0,one = 1.0;
998: comm = eps->comm;
1000: MPI_Comm_size(comm,&size);
1002: VecDuplicate(eps->vec_initial,&in);
1003: VecDuplicate(eps->vec_initial,&out);
1004: VecGetSize(in,&M);
1005: VecGetLocalSize(in,&m);
1006: VecGetOwnershipRange(in,&start,&end);
1007: PetscMalloc((m+1)*sizeof(int),&rows);
1008: for (i=0; i<m; i++) {rows[i] = start + i;}
1010: if (size == 1) {
1011: MatCreateSeqDense(comm,M,M,PETSC_NULL,mat);
1012: } else {
1013: MatCreateMPIAIJ(comm,m,m,M,M,0,0,0,0,mat);
1014: }
1015:
1016: for (i=0; i<M; i++) {
1018: VecSet(&zero,in);
1019: VecSetValues(in,1,&i,&one,INSERT_VALUES);
1020: VecAssemblyBegin(in);
1021: VecAssemblyEnd(in);
1023: STApply(eps->OP,in,out);
1024:
1025: VecGetArray(out,&array);
1026: MatSetValues(*mat,m,rows,1,&i,array,INSERT_VALUES);
1027: VecRestoreArray(out,&array);
1029: }
1030: PetscFree(rows);
1031: VecDestroy(in);
1032: VecDestroy(out);
1033: MatAssemblyBegin(*mat,MAT_FINAL_ASSEMBLY);
1034: MatAssemblyEnd(*mat,MAT_FINAL_ASSEMBLY);
1035: return(0);
1036: }
1040: /*@
1041: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
1043: Collective on EPS and Mat
1045: Input Parameters:
1046: + eps - the eigenproblem solver context
1047: . A - the matrix associated with the eigensystem
1048: - B - the second matrix in the case of generalized eigenproblems
1050: Notes:
1051: To specify a standard eigenproblem, use PETSC_NULL for parameter B.
1053: Level: beginner
1055: .seealso: EPSSolve(), EPSGetST(), STGetOperators()
1056: @*/
1057: int EPSSetOperators(EPS eps,Mat A,Mat B)
1058: {
1065: STSetOperators(eps->OP,A,B);
1066: eps->setupcalled = 0; /* so that next solve call will call setup */
1068: /* The following call is done in order to check the consistency of the
1069: problem type with the specified matrices */
1070: if (eps->problem_type) {
1071: EPSSetProblemType(eps,eps->problem_type);
1072: }
1073: return(0);
1074: }
1078: /*@
1079: EPSComputeResidualNorm - Computes the residual norm associated with
1080: the i-th converged approximate eigenpair.
1082: Collective on EPS
1084: Input Parameter:
1085: . eps - the eigensolver context
1086: . i - the solution index
1088: Output Parameter:
1089: . norm - the residual norm, computed as ||Ax-kBx|| where k is the
1090: eigenvalue and x is the eigenvector.
1091: If k=0 then the residual norm is computed as ||Ax||.
1093: Notes:
1094: The index i should be a value between 0 and nconv (see EPSGetConverged()).
1095: Eigenpairs are indexed according to the ordering criterion established
1096: with EPSSetWhichEigenpairs().
1098: Level: beginner
1100: .seealso: EPSSolve(), EPSGetConverged(), EPSSetWhichEigenpairs()
1101: @*/
1102: int EPSComputeResidualNorm(EPS eps, int i, PetscReal *norm)
1103: {
1104: Vec u, v, w, xr, xi;
1105: Mat A, B;
1106: int ierr;
1107: PetscScalar alpha, kr, ki;
1108: PetscReal ni, nr;
1109:
1112: if (eps->dropvectors || !eps->V) { SETERRQ(1, "Eigenvectors are not available"); }
1113: STGetOperators(eps->OP,&A,&B);
1114: VecDuplicate(eps->vec_initial,&u);
1115: VecDuplicate(eps->vec_initial,&v);
1116: VecDuplicate(eps->vec_initial,&w);
1117: VecDuplicate(eps->vec_initial,&xr);
1118: VecDuplicate(eps->vec_initial,&xi);
1119: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
1121: #ifndef PETSC_USE_COMPLEX
1122: if (ki == 0 ||
1123: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1124: #endif
1125: MatMult( A, xr, u ); /* u=A*x */
1126: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1127: if (eps->isgeneralized) { MatMult( B, xr, w ); }
1128: else { VecCopy( xr, w ); } /* w=B*x */
1129: alpha = -kr;
1130: VecAXPY( &alpha, w, u ); /* u=A*x-k*B*x */
1131: }
1132: VecNorm( u, NORM_2, norm);
1133: #ifndef PETSC_USE_COMPLEX
1134: } else {
1135: MatMult( A, xr, u ); /* u=A*xr */
1136: if (eps->isgeneralized) { MatMult( B, xr, v ); }
1137: else { VecCopy( xr, v ); } /* v=B*xr */
1138: alpha = -kr;
1139: VecAXPY( &alpha, v, u ); /* u=A*xr-kr*B*xr */
1140: if (eps->isgeneralized) { MatMult( B, xi, w ); }
1141: else { VecCopy( xi, w ); } /* w=B*xi */
1142: alpha = ki;
1143: VecAXPY( &alpha, w, u ); /* u=A*xr-kr*B*xr+ki*B*xi */
1144: VecNorm( u, NORM_2, &nr );
1145: MatMult( A, xi, u ); /* u=A*xi */
1146: alpha = -kr;
1147: VecAXPY( &alpha, w, u ); /* u=A*xi-kr*B*xi */
1148: alpha = -ki;
1149: VecAXPY( &alpha, v, u ); /* u=A*xi-kr*B*xi-ki*B*xr */
1150: VecNorm( u, NORM_2, &ni );
1151: *norm = LAlapy2_( &nr, &ni );
1152: }
1153: #endif
1155: VecDestroy(w);
1156: VecDestroy(v);
1157: VecDestroy(u);
1158: VecDestroy(xr);
1159: VecDestroy(xi);
1160: return(0);
1161: }
1165: /*@
1166: EPSComputeRelativeError - Computes the actual relative error associated
1167: with the i-th converged approximate eigenpair.
1169: Collective on EPS
1171: Input Parameter:
1172: . eps - the eigensolver context
1173: . i - the solution index
1175: Output Parameter:
1176: . error - the relative error, computed as ||Ax-kBx||/||kx|| where k is the
1177: eigenvalue and x is the eigenvector.
1178: If k=0 the relative error is computed as ||Ax||/||x||.
1180: Level: beginner
1182: .seealso: EPSSolve(), EPSComputeResidualNorm()
1183: @*/
1184: int EPSComputeRelativeError(EPS eps, int i, PetscReal *error)
1185: {
1186: Vec xr, xi, u;
1187: int ierr;
1188: PetscScalar kr, ki, alpha;
1189: PetscReal norm, er, ei;
1190:
1193: EPSComputeResidualNorm(eps,i,&norm);
1194: VecDuplicate(eps->vec_initial,&xr);
1195: VecDuplicate(eps->vec_initial,&xi);
1196: EPSGetEigenpair(eps,i,&kr,&ki,xr,xi);
1198: #ifndef PETSC_USE_COMPLEX
1199: if (ki == 0 ||
1200: PetscAbsScalar(ki) < PetscAbsScalar(kr*PETSC_MACHINE_EPSILON)) {
1201: #endif
1202: if (PetscAbsScalar(kr) > PETSC_MACHINE_EPSILON) {
1203: VecScale(&kr, xr);
1204: }
1205: VecNorm(xr, NORM_2, &er);
1206: *error = norm / er;
1207: #ifndef PETSC_USE_COMPLEX
1208: } else {
1209: VecDuplicate(xi, &u);
1210: VecCopy(xi, u);
1211: alpha = -ki;
1212: VecAXPBY(&kr, &alpha, xr, u);
1213: VecNorm(u, NORM_2, &er);
1214: VecAXPBY(&kr, &ki, xr, xi);
1215: VecNorm(xi, NORM_2, &ei);
1216: VecDestroy(u);
1217: *error = norm / LAlapy2_(&er, &ei);
1218: }
1219: #endif
1220:
1221: VecDestroy(xr);
1222: VecDestroy(xi);
1223: return(0);
1224: }
1228: /*@
1229: EPSSetProblemType - Specifies the type of the eigenvalue problem.
1231: Collective on EPS
1233: Input Parameters:
1234: + eps - the eigensolver context
1235: - type - a known type of eigenvalue problem
1237: Options Database Keys:
1238: + -eps_hermitian - Hermitian eigenvalue problem (default)
1239: . -eps_gen_hermitian - generalized Hermitian eigenvalue problem
1240: . -eps_non_hermitian - non-Hermitian eigenvalue problem
1241: - -eps_gen_non_hermitian - generalized non-Hermitian eigenvalue problem
1242:
1243: Note:
1244: Normally, the user need not set the EPS type, since it can be determined from
1245: the information given in the EPSSetOperators call. This routine is reserved
1246: for special cases such as when a nonsymmetric solver wants to be
1247: used in a symmetric problem.
1249: Level: advanced
1251: .seealso: EPSSetOperators(), EPSSetType(), EPSType
1252: @*/
1253: int EPSSetProblemType(EPS eps,EPSProblemType type)
1254: {
1255: int n,m,ierr;
1256: Mat A,B;
1257: PetscTruth Ah,Bh,inconsistent=PETSC_FALSE;
1262: if (type!=EPS_HEP && type!=EPS_GHEP && type!=EPS_NHEP && type!=EPS_GNHEP ) { SETERRQ(PETSC_ERR_ARG_WRONG,"Unknown eigenvalue problem type"); }
1264: STGetOperators(eps->OP,&A,&B);
1265: if (!A) { SETERRQ(1,"Must call EPSSetOperators() first"); }
1267: /* Check for square matrices */
1268: MatGetSize(A,&m,&n);
1269: if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
1270: if (B) {
1271: MatGetSize(B,&m,&n);
1272: if (m!=n) { SETERRQ(1,"B is a non-square matrix"); }
1273: }
1275: eps->problem_type = type;
1277: SlepcIsHermitian(A,&Ah);
1278: if (B) { SlepcIsHermitian(B,&Bh); }
1280: if (!B) {
1281: eps->isgeneralized = PETSC_FALSE;
1282: if (Ah) eps->ishermitian = PETSC_TRUE;
1283: else eps->ishermitian = PETSC_FALSE;
1284: }
1285: else {
1286: eps->isgeneralized = PETSC_TRUE;
1287: if (Ah && Bh) eps->ishermitian = PETSC_TRUE;
1288: else eps->ishermitian = PETSC_FALSE;
1289: }
1290:
1291: switch (type) {
1292: case EPS_HEP:
1293: if (eps->isgeneralized || !eps->ishermitian) inconsistent=PETSC_TRUE;
1294: eps->ishermitian = PETSC_TRUE;
1295: break;
1296: case EPS_GHEP:
1297: /* Note that here we do not consider the case in which A and B are
1298: non-hermitian but there exists a linear combination of them which is */
1299: if (!eps->isgeneralized || !eps->ishermitian) inconsistent=PETSC_TRUE;
1300: break;
1301: case EPS_NHEP:
1302: if (eps->isgeneralized) inconsistent=PETSC_TRUE;
1303: eps->ishermitian = PETSC_FALSE;
1304: break;
1305: case EPS_GNHEP:
1306: /* If matrix B is not given then an error is issued. An alternative
1307: would be to generate an identity matrix. Also in EPS_GHEP above */
1308: if (!eps->isgeneralized) inconsistent=PETSC_TRUE;
1309: eps->ishermitian = PETSC_FALSE;
1310: break;
1311: }
1312: if (inconsistent) { SETERRQ(0,"Warning: Inconsistent EPS state"); }
1314: return(0);
1315: }
1319: /*@
1320: EPSGetProblemType - Gets the problem type from the EPS object.
1322: Not Collective
1324: Input Parameter:
1325: . eps - the eigensolver context
1327: Output Parameter:
1328: . type - name of EPS problem type
1330: Level: intermediate
1332: .seealso: EPSSetProblemType()
1333: @*/
1334: int EPSGetProblemType(EPS eps,EPSProblemType *type)
1335: {
1338: *type = eps->problem_type;
1339: return(0);
1340: }
1344: /*@
1345: EPSIsGeneralized - Ask if the EPS object corresponds to a generalized
1346: eigenvalue problem.
1348: Not collective
1350: Input Parameter:
1351: . eps - the eigenproblem solver context
1353: Output Parameter:
1354: . is - the answer
1356: Level: intermediate
1358: @*/
1359: int EPSIsGeneralized(EPS eps,PetscTruth* is)
1360: {
1361: int ierr;
1362: Mat B;
1366: STGetOperators(eps->OP,PETSC_NULL,&B);
1367: if( B ) *is = PETSC_TRUE;
1368: else *is = PETSC_FALSE;
1369: if( eps->setupcalled ) {
1370: if( eps->isgeneralized != *is ) {
1371: SETERRQ(0,"Warning: Inconsistent EPS state");
1372: }
1373: }
1374: return(0);
1375: }
1379: /*@
1380: EPSIsHermitian - Ask if the EPS object corresponds to a Hermitian
1381: eigenvalue problem.
1383: Not collective
1385: Input Parameter:
1386: . eps - the eigenproblem solver context
1388: Output Parameter:
1389: . is - the answer
1391: Level: intermediate
1393: @*/
1394: int EPSIsHermitian(EPS eps,PetscTruth* is)
1395: {
1398: if( eps->ishermitian ) *is = PETSC_TRUE;
1399: else *is = PETSC_FALSE;
1400: return(0);
1401: }
1405: /*@
1406: EPSReverseProjection - Compute the operation V=V*S, where the columns of
1407: V are m of the basis vectors of the EPS object and S is an mxm dense
1408: matrix.
1410: Collective on EPS
1412: Input Parameter:
1413: + eps - the eigenproblem solver context
1414: . k - starting column
1415: . m - dimension of matrix S
1416: - S - pointer to the values of matrix S
1418: Level: developer
1420: Note:
1421: Matrix S is overwritten.
1423: @*/
1424: int EPSReverseProjection(EPS eps,int k,int m,PetscScalar *S)
1425: {
1426: int i,j,n,ierr,lwork;
1427: PetscScalar *tau,*work,*pV;
1428:
1431: VecGetLocalSize(eps->vec_initial,&n);
1432: lwork = n;
1433: PetscMalloc(m*sizeof(PetscScalar),&tau);
1434: PetscMalloc(lwork*sizeof(PetscScalar),&work);
1436: /* compute the LQ factorization L Q = S */
1437: LAgelqf_(&m,&m,S,&m,tau,work,&lwork,&ierr);
1439: /* triangular post-multiplication, V = V L */
1440: for (i=k;i<k+m;i++) {
1441: VecScale(S+(i-k)+m*(i-k),eps->V[i]);
1442: for (j=i+1;j<k+m;j++) {
1443: VecAXPY(S+(j-k)+m*(i-k),eps->V[j],eps->V[i]);
1444: }
1445: }
1447: /* orthogonal post-multiplication, V = V Q */
1448: VecGetArray(eps->V[k],&pV);
1449: LAormlq_("R","N",&n,&m,&m,S,&m,tau,pV,&n,work,&lwork,&ierr,1,1);
1450: VecRestoreArray(eps->V[k],&pV);
1452: PetscFree(tau);
1453: PetscFree(work);
1455: return(0);
1456: }
1461: /*@
1462: EPSSwapEigenpairs - Swaps all the information internal to the EPS object
1463: corresponding to eigenpairs which occupy the i-th and j-th positions.
1465: Collective on EPS
1467: Input Parameter:
1468: + eps - the eigenproblem solver context
1469: . i - first index
1470: - j - second index
1472: Level: developer
1474: @*/
1475: int EPSSwapEigenpairs(EPS eps,int i,int j)
1476: {
1477: int ierr;
1478: PetscScalar tscalar;
1479: PetscReal treal;
1480:
1482: if (i!=j) {
1483: VecSwap(eps->V[i],eps->V[j]);
1484: tscalar = eps->eigr[i];
1485: eps->eigr[i] = eps->eigr[j];
1486: eps->eigr[j] = tscalar;
1487: tscalar = eps->eigi[i];
1488: eps->eigi[i] = eps->eigi[j];
1489: eps->eigi[j] = tscalar;
1490: treal = eps->errest[i];
1491: eps->errest[i] = eps->errest[j];
1492: eps->errest[j] = treal;
1493: }
1494: return(0);
1495: }
1499: int EPSBackTransform_Default(EPS eps)
1500: {
1501: ST st;
1502: int ierr,i;
1506: EPSGetST(eps,&st);
1507: for (i=0;i<eps->nconv;i++) {
1508: STBackTransform(st,&eps->eigr[i],&eps->eigi[i]);
1509: }
1511: return(0);
1512: }