1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain
6: This file is part of SLEPc.
7: SLEPc is distributed under a 2-clause BSD license (see LICENSE).
8: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
9: */
10: /*
11: EPS routines related to problem setup
12: */
14: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
16: /*
17: Let the solver choose the ST type that should be used by default,
18: otherwise set it to SHIFT.
19: This is called at EPSSetFromOptions (before STSetFromOptions)
20: and also at EPSSetUp (in case EPSSetFromOptions was not called).
21: */
22: PetscErrorCode EPSSetDefaultST(EPS eps) 23: {
27: if (eps->ops->setdefaultst) { (*eps->ops->setdefaultst)(eps); }
28: if (!((PetscObject)eps->st)->type_name) {
29: STSetType(eps->st,STSHIFT);
30: }
31: return(0);
32: }
34: /*
35: This is done by preconditioned eigensolvers that use the PC only.
36: It sets STPRECOND with KSPPREONLY.
37: */
38: PetscErrorCode EPSSetDefaultST_Precond(EPS eps) 39: {
41: KSP ksp;
44: if (!((PetscObject)eps->st)->type_name) {
45: STSetType(eps->st,STPRECOND);
46: STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
47: }
48: STGetKSP(eps->st,&ksp);
49: if (!((PetscObject)ksp)->type_name) {
50: KSPSetType(ksp,KSPPREONLY);
51: }
52: return(0);
53: }
55: /*
56: This is done by preconditioned eigensolvers that can also use the KSP.
57: It sets STPRECOND with the default KSP (GMRES).
58: */
59: PetscErrorCode EPSSetDefaultST_GMRES(EPS eps) 60: {
62: KSP ksp;
65: if (!((PetscObject)eps->st)->type_name) {
66: STSetType(eps->st,STPRECOND);
67: STPrecondSetKSPHasMat(eps->st,PETSC_TRUE);
68: }
69: STGetKSP(eps->st,&ksp);
70: KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,5);
71: return(0);
72: }
74: /*@
75: EPSSetUp - Sets up all the internal data structures necessary for the
76: execution of the eigensolver. Then calls STSetUp() for any set-up
77: operations associated to the ST object.
79: Collective on EPS 81: Input Parameter:
82: . eps - eigenproblem solver context
84: Notes:
85: This function need not be called explicitly in most cases, since EPSSolve()
86: calls it. It can be useful when one wants to measure the set-up time
87: separately from the solve time.
89: Level: developer
91: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp(), EPSSetInitialSpace()
92: @*/
93: PetscErrorCode EPSSetUp(EPS eps) 94: {
96: Mat A,B;
97: SlepcSC sc;
98: PetscInt k,nmat;
99: PetscBool flg,istrivial,precond;
100: #if defined(PETSC_USE_COMPLEX)
101: PetscScalar sigma;
102: #endif
106: if (eps->state) return(0);
107: PetscLogEventBegin(EPS_SetUp,eps,0,0,0);
109: /* reset the convergence flag from the previous solves */
110: eps->reason = EPS_CONVERGED_ITERATING;
112: /* Set default solver type (EPSSetFromOptions was not called) */
113: if (!((PetscObject)eps)->type_name) {
114: EPSSetType(eps,EPSKRYLOVSCHUR);
115: }
116: if (!eps->st) { EPSGetST(eps,&eps->st); }
117: EPSSetDefaultST(eps);
119: PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&precond);
120: if (eps->categ==EPS_CATEGORY_PRECOND && !precond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires ST=PRECOND");
121: if (eps->categ!=EPS_CATEGORY_PRECOND && precond) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"STPRECOND is intended for preconditioned eigensolvers only");
123: STSetTransform(eps->st,PETSC_TRUE);
124: if (eps->useds && !eps->ds) { EPSGetDS(eps,&eps->ds); }
125: if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
126: if (!((PetscObject)eps->rg)->type_name) {
127: RGSetType(eps->rg,RGINTERVAL);
128: }
130: /* Set problem dimensions */
131: STGetNumMatrices(eps->st,&nmat);
132: if (!nmat) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
133: STMatGetSize(eps->st,&eps->n,NULL);
134: STMatGetLocalSize(eps->st,&eps->nloc,NULL);
136: /* Set default problem type */
137: if (!eps->problem_type) {
138: if (nmat==1) {
139: EPSSetProblemType(eps,EPS_NHEP);
140: } else {
141: EPSSetProblemType(eps,EPS_GNHEP);
142: }
143: } else if (nmat==1 && eps->isgeneralized) {
144: PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n");
145: eps->isgeneralized = PETSC_FALSE;
146: eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
147: } else if (nmat>1 && !eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state");
149: if (eps->nev > eps->n) eps->nev = eps->n;
150: if (eps->ncv > eps->n) eps->ncv = eps->n;
152: /* initialization of matrix norms */
153: if (eps->conv==EPS_CONV_NORM) {
154: if (!eps->nrma) {
155: STGetMatrix(eps->st,0,&A);
156: MatNorm(A,NORM_INFINITY,&eps->nrma);
157: }
158: if (nmat>1 && !eps->nrmb) {
159: STGetMatrix(eps->st,1,&B);
160: MatNorm(B,NORM_INFINITY,&eps->nrmb);
161: }
162: }
164: /* call specific solver setup */
165: (*eps->ops->setup)(eps);
167: /* if purification is set, check that it really makes sense */
168: if (eps->purify) {
169: if (eps->categ==EPS_CATEGORY_PRECOND || eps->categ==EPS_CATEGORY_CONTOUR) eps->purify = PETSC_FALSE;
170: else {
171: if (!eps->isgeneralized) eps->purify = PETSC_FALSE;
172: else if (!eps->ishermitian && !eps->ispositive) eps->purify = PETSC_FALSE;
173: else {
174: PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg);
175: if (flg) eps->purify = PETSC_FALSE;
176: }
177: }
178: }
180: /* check extraction */
181: PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STPRECOND,STSHIFT,"");
182: if (!flg && eps->extraction && eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot use a spectral transformation combined with harmonic extraction");
184: /* set tolerance if not yet set */
185: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
187: /* fill sorting criterion context */
188: switch (eps->which) {
189: case EPS_LARGEST_MAGNITUDE:
190: eps->sc->comparison = SlepcCompareLargestMagnitude;
191: eps->sc->comparisonctx = NULL;
192: break;
193: case EPS_SMALLEST_MAGNITUDE:
194: eps->sc->comparison = SlepcCompareSmallestMagnitude;
195: eps->sc->comparisonctx = NULL;
196: break;
197: case EPS_LARGEST_REAL:
198: eps->sc->comparison = SlepcCompareLargestReal;
199: eps->sc->comparisonctx = NULL;
200: break;
201: case EPS_SMALLEST_REAL:
202: eps->sc->comparison = SlepcCompareSmallestReal;
203: eps->sc->comparisonctx = NULL;
204: break;
205: case EPS_LARGEST_IMAGINARY:
206: eps->sc->comparison = SlepcCompareLargestImaginary;
207: eps->sc->comparisonctx = NULL;
208: break;
209: case EPS_SMALLEST_IMAGINARY:
210: eps->sc->comparison = SlepcCompareSmallestImaginary;
211: eps->sc->comparisonctx = NULL;
212: break;
213: case EPS_TARGET_MAGNITUDE:
214: eps->sc->comparison = SlepcCompareTargetMagnitude;
215: eps->sc->comparisonctx = &eps->target;
216: break;
217: case EPS_TARGET_REAL:
218: eps->sc->comparison = SlepcCompareTargetReal;
219: eps->sc->comparisonctx = &eps->target;
220: break;
221: case EPS_TARGET_IMAGINARY:
222: #if defined(PETSC_USE_COMPLEX)
223: eps->sc->comparison = SlepcCompareTargetImaginary;
224: eps->sc->comparisonctx = &eps->target;
225: #endif
226: break;
227: case EPS_ALL:
228: eps->sc->comparison = SlepcCompareSmallestReal;
229: eps->sc->comparisonctx = NULL;
230: break;
231: case EPS_WHICH_USER:
232: break;
233: }
234: eps->sc->map = NULL;
235: eps->sc->mapobj = NULL;
237: if (eps->useds) {
238: /* fill sorting criterion for DS */
239: DSGetSlepcSC(eps->ds,&sc);
240: RGIsTrivial(eps->rg,&istrivial);
241: if (eps->which!=EPS_ALL) {
242: sc->rg = istrivial? NULL: eps->rg;
243: sc->comparison = eps->sc->comparison;
244: sc->comparisonctx = eps->sc->comparisonctx;
245: sc->map = SlepcMap_ST;
246: sc->mapobj = (PetscObject)eps->st;
247: }
248: }
250: /* Build balancing matrix if required */
251: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
252: if (!eps->D) {
253: BVCreateVec(eps->V,&eps->D);
254: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->D);
255: } else {
256: VecSet(eps->D,1.0);
257: }
258: EPSBuildBalance_Krylov(eps);
259: STSetBalanceMatrix(eps->st,eps->D);
260: }
262: /* Setup ST */
263: STSetUp(eps->st);
265: #if defined(PETSC_USE_COMPLEX)
266: STGetShift(eps->st,&sigma);
267: if (eps->ishermitian && PetscImaginaryPart(sigma) != 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
268: #endif
269: PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg);
270: if (flg && eps->problem_type == EPS_PGNHEP) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
271: PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STSINVERT,STCAYLEY,"");
272: if (flg && (eps->categ==EPS_CATEGORY_KRYLOV || eps->categ==EPS_CATEGORY_OTHER) && (eps->which!=EPS_TARGET_MAGNITUDE && eps->which!=EPS_TARGET_REAL && eps->which!=EPS_TARGET_IMAGINARY && eps->which!=EPS_ALL && eps->which!=EPS_WHICH_USER)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"Shift-and-invert requires a target 'which' (see EPSSetWhichEigenpairs), for instance -st_type sinvert -eps_target 0 -eps_target_magnitude");
274: /* process deflation and initial vectors */
275: if (eps->nds<0) {
276: k = -eps->nds;
277: BVInsertConstraints(eps->V,&k,eps->defl);
278: SlepcBasisDestroy_Private(&eps->nds,&eps->defl);
279: eps->nds = k;
280: STCheckNullSpace(eps->st,eps->V);
281: }
282: if (eps->nini<0) {
283: k = -eps->nini;
284: if (k>eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The number of initial vectors is larger than ncv");
285: BVInsertVecs(eps->V,0,&k,eps->IS,PETSC_TRUE);
286: SlepcBasisDestroy_Private(&eps->nini,&eps->IS);
287: eps->nini = k;
288: }
290: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
291: eps->state = EPS_STATE_SETUP;
292: return(0);
293: }
295: /*@C
296: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
298: Collective on EPS and Mat
300: Input Parameters:
301: + eps - the eigenproblem solver context
302: . A - the matrix associated with the eigensystem
303: - B - the second matrix in the case of generalized eigenproblems
305: Notes:
306: To specify a standard eigenproblem, use NULL for parameter B.
308: It must be called before EPSSetUp(). If it is called again after EPSSetUp() then
309: the EPS object is reset.
311: Level: beginner
313: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetMatrix()
314: @*/
315: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)316: {
318: PetscInt m,n,m0,nmat;
319: Mat mat[2];
328: /* Check for square matrices */
329: MatGetSize(A,&m,&n);
330: if (m!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix");
331: if (B) {
332: MatGetSize(B,&m0,&n);
333: if (m0!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix");
334: if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match");
335: }
336: if (eps->state && n!=eps->n) { EPSReset(eps); }
337: eps->nrma = 0.0;
338: eps->nrmb = 0.0;
339: if (!eps->st) { EPSGetST(eps,&eps->st); }
340: mat[0] = A;
341: if (B) {
342: mat[1] = B;
343: nmat = 2;
344: } else nmat = 1;
345: STSetMatrices(eps->st,nmat,mat);
346: eps->state = EPS_STATE_INITIAL;
347: return(0);
348: }
350: /*@
351: EPSGetOperators - Gets the matrices associated with the eigensystem.
353: Collective on EPS and Mat
355: Input Parameter:
356: . eps - the EPS context
358: Output Parameters:
359: + A - the matrix associated with the eigensystem
360: - B - the second matrix in the case of generalized eigenproblems
362: Level: intermediate
364: .seealso: EPSSolve(), EPSGetST(), STGetMatrix(), STSetMatrices()
365: @*/
366: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)367: {
369: ST st;
370: PetscInt k;
374: EPSGetST(eps,&st);
375: if (A) { STGetMatrix(st,0,A); }
376: if (B) {
377: STGetNumMatrices(st,&k);
378: if (k==1) B = NULL;
379: else {
380: STGetMatrix(st,1,B);
381: }
382: }
383: return(0);
384: }
386: /*@
387: EPSSetDeflationSpace - Specify a basis of vectors that constitute the deflation
388: space.
390: Collective on EPS and Vec
392: Input Parameter:
393: + eps - the eigenproblem solver context
394: . n - number of vectors
395: - v - set of basis vectors of the deflation space
397: Notes:
398: When a deflation space is given, the eigensolver seeks the eigensolution
399: in the restriction of the problem to the orthogonal complement of this
400: space. This can be used for instance in the case that an invariant
401: subspace is known beforehand (such as the nullspace of the matrix).
403: These vectors do not persist from one EPSSolve() call to the other, so the
404: deflation space should be set every time.
406: The vectors do not need to be mutually orthonormal, since they are explicitly
407: orthonormalized internally.
409: Level: intermediate
411: .seealso: EPSSetInitialSpace()
412: @*/
413: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *v)414: {
420: if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
421: SlepcBasisReference_Private(n,v,&eps->nds,&eps->defl);
422: if (n>0) eps->state = EPS_STATE_INITIAL;
423: return(0);
424: }
426: /*@C
427: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
428: space, that is, the subspace from which the solver starts to iterate.
430: Collective on EPS and Vec
432: Input Parameter:
433: + eps - the eigenproblem solver context
434: . n - number of vectors
435: - is - set of basis vectors of the initial space
437: Notes:
438: Some solvers start to iterate on a single vector (initial vector). In that case,
439: the other vectors are ignored.
441: These vectors do not persist from one EPSSolve() call to the other, so the
442: initial space should be set every time.
444: The vectors do not need to be mutually orthonormal, since they are explicitly
445: orthonormalized internally.
447: Common usage of this function is when the user can provide a rough approximation
448: of the wanted eigenspace. Then, convergence may be faster.
450: Level: intermediate
452: .seealso: EPSSetDeflationSpace()
453: @*/
454: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)455: {
461: if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
462: SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS);
463: if (n>0) eps->state = EPS_STATE_INITIAL;
464: return(0);
465: }
467: /*
468: EPSSetDimensions_Default - Set reasonable values for ncv, mpd if not set
469: by the user. This is called at setup.
470: */
471: PetscErrorCode EPSSetDimensions_Default(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)472: {
474: PetscBool krylov;
477: if (*ncv) { /* ncv set */
478: PetscObjectTypeCompareAny((PetscObject)eps,&krylov,EPSKRYLOVSCHUR,EPSARNOLDI,EPSLANCZOS,"");
479: if (krylov) {
480: if (*ncv<nev+1 && !(*ncv==nev && *ncv==eps->n)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev+1");
481: } else {
482: if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
483: }
484: } else if (*mpd) { /* mpd set */
485: *ncv = PetscMin(eps->n,nev+(*mpd));
486: } else { /* neither set: defaults depend on nev being small or large */
487: if (nev<500) *ncv = PetscMin(eps->n,PetscMax(2*nev,nev+15));
488: else {
489: *mpd = 500;
490: *ncv = PetscMin(eps->n,nev+(*mpd));
491: }
492: }
493: if (!*mpd) *mpd = *ncv;
494: return(0);
495: }
497: /*@
498: EPSAllocateSolution - Allocate memory storage for common variables such
499: as eigenvalues and eigenvectors.
501: Collective on EPS503: Input Parameters:
504: + eps - eigensolver context
505: - extra - number of additional positions, used for methods that require a
506: working basis slightly larger than ncv
508: Developers Note:
509: This is PETSC_EXTERN because it may be required by user plugin EPS510: implementations.
512: Level: developer
513: @*/
514: PetscErrorCode EPSAllocateSolution(EPS eps,PetscInt extra)515: {
517: PetscInt oldsize,newc,requested;
518: PetscLogDouble cnt;
519: Vec t;
522: requested = eps->ncv + extra;
524: /* oldsize is zero if this is the first time setup is called */
525: BVGetSizes(eps->V,NULL,NULL,&oldsize);
526: newc = PetscMax(0,requested-oldsize);
528: /* allocate space for eigenvalues and friends */
529: if (requested != oldsize || !eps->eigr) {
530: PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm);
531: PetscMalloc4(requested,&eps->eigr,requested,&eps->eigi,requested,&eps->errest,requested,&eps->perm);
532: cnt = 2*newc*sizeof(PetscScalar) + 2*newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
533: PetscLogObjectMemory((PetscObject)eps,cnt);
534: }
536: /* workspace for the case of arbitrary selection */
537: if (eps->arbitrary) {
538: if (eps->rr) {
539: PetscFree2(eps->rr,eps->ri);
540: }
541: PetscMalloc2(requested,&eps->rr,requested,&eps->ri);
542: PetscLogObjectMemory((PetscObject)eps,2*newc*sizeof(PetscScalar));
543: }
545: /* allocate V */
546: if (!eps->V) { EPSGetBV(eps,&eps->V); }
547: if (!oldsize) {
548: if (!((PetscObject)(eps->V))->type_name) {
549: BVSetType(eps->V,BVSVEC);
550: }
551: STMatCreateVecsEmpty(eps->st,&t,NULL);
552: BVSetSizesFromVec(eps->V,t,requested);
553: VecDestroy(&t);
554: } else {
555: BVResize(eps->V,requested,PETSC_FALSE);
556: }
557: return(0);
558: }