1: /*
2: EPS routines related to problem setup.
4: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
5: SLEPc - Scalable Library for Eigenvalue Problem Computations
6: Copyright (c) 2002-2015, Universitat Politecnica de Valencia, Spain
8: This file is part of SLEPc.
10: SLEPc is free software: you can redistribute it and/or modify it under the
11: terms of version 3 of the GNU Lesser General Public License as published by
12: the Free Software Foundation.
14: SLEPc is distributed in the hope that it will be useful, but WITHOUT ANY
15: WARRANTY; without even the implied warranty of MERCHANTABILITY or FITNESS
16: FOR A PARTICULAR PURPOSE. See the GNU Lesser General Public License for
17: more details.
19: You should have received a copy of the GNU Lesser General Public License
20: along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
21: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
22: */
24: #include <slepc/private/epsimpl.h> /*I "slepceps.h" I*/
28: /*@
29: EPSSetUp - Sets up all the internal data structures necessary for the
30: execution of the eigensolver. Then calls STSetUp() for any set-up
31: operations associated to the ST object.
33: Collective on EPS 35: Input Parameter:
36: . eps - eigenproblem solver context
38: Notes:
39: This function need not be called explicitly in most cases, since EPSSolve()
40: calls it. It can be useful when one wants to measure the set-up time
41: separately from the solve time.
43: Level: developer
45: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp(), EPSSetInitialSpace()
46: @*/
47: PetscErrorCode EPSSetUp(EPS eps) 48: {
50: Mat A,B;
51: SlepcSC sc;
52: PetscInt k,nmat;
53: PetscBool flg,istrivial;
54: #if defined(PETSC_USE_COMPLEX)
55: PetscScalar sigma;
56: #endif
60: if (eps->state) return(0);
61: PetscLogEventBegin(EPS_SetUp,eps,0,0,0);
63: /* reset the convergence flag from the previous solves */
64: eps->reason = EPS_CONVERGED_ITERATING;
66: /* Set default solver type (EPSSetFromOptions was not called) */
67: if (!((PetscObject)eps)->type_name) {
68: EPSSetType(eps,EPSKRYLOVSCHUR);
69: }
70: if (!eps->st) { EPSGetST(eps,&eps->st); }
71: if (!((PetscObject)eps->st)->type_name) {
72: PetscObjectTypeCompareAny((PetscObject)eps,&flg,EPSGD,EPSJD,EPSRQCG,EPSBLOPEX,EPSPRIMME,"");
73: STSetType(eps->st,flg?STPRECOND:STSHIFT);
74: }
75: STSetTransform(eps->st,PETSC_TRUE);
76: if (!eps->ds) { EPSGetDS(eps,&eps->ds); }
77: DSReset(eps->ds);
78: if (!eps->rg) { EPSGetRG(eps,&eps->rg); }
79: if (!((PetscObject)eps->rg)->type_name) {
80: RGSetType(eps->rg,RGINTERVAL);
81: }
82: if (!((PetscObject)eps->rand)->type_name) {
83: PetscRandomSetFromOptions(eps->rand);
84: }
86: /* Set problem dimensions */
87: STGetNumMatrices(eps->st,&nmat);
88: if (!nmat) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
89: STMatGetSize(eps->st,&eps->n,NULL);
90: STMatGetLocalSize(eps->st,&eps->nloc,NULL);
92: /* Set default problem type */
93: if (!eps->problem_type) {
94: if (nmat==1) {
95: EPSSetProblemType(eps,EPS_NHEP);
96: } else {
97: EPSSetProblemType(eps,EPS_GNHEP);
98: }
99: } else if (nmat==1 && eps->isgeneralized) {
100: PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n");
101: eps->isgeneralized = PETSC_FALSE;
102: eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
103: } else if (nmat>1 && !eps->isgeneralized) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state");
105: if (eps->nev > eps->n) eps->nev = eps->n;
106: if (eps->ncv > eps->n) eps->ncv = eps->n;
108: /* initialization of matrix norms */
109: if (eps->conv==EPS_CONV_NORM) {
110: if (!eps->nrma) {
111: STGetOperators(eps->st,0,&A);
112: MatNorm(A,NORM_INFINITY,&eps->nrma);
113: }
114: if (nmat>1 && !eps->nrmb) {
115: STGetOperators(eps->st,1,&B);
116: MatNorm(B,NORM_INFINITY,&eps->nrmb);
117: }
118: }
120: /* call specific solver setup */
121: (*eps->ops->setup)(eps);
123: /* check extraction */
124: PetscObjectTypeCompareAny((PetscObject)eps->st,&flg,STPRECOND,STSHIFT,"");
125: if (!flg && eps->extraction && eps->extraction!=EPS_RITZ) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot use a spectral transformation combined with harmonic extraction");
127: /* set tolerance if not yet set */
128: if (eps->tol==PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
130: /* fill sorting criterion context */
131: switch (eps->which) {
132: case EPS_LARGEST_MAGNITUDE:
133: eps->sc->comparison = SlepcCompareLargestMagnitude;
134: eps->sc->comparisonctx = NULL;
135: break;
136: case EPS_SMALLEST_MAGNITUDE:
137: eps->sc->comparison = SlepcCompareSmallestMagnitude;
138: eps->sc->comparisonctx = NULL;
139: break;
140: case EPS_LARGEST_REAL:
141: eps->sc->comparison = SlepcCompareLargestReal;
142: eps->sc->comparisonctx = NULL;
143: break;
144: case EPS_SMALLEST_REAL:
145: eps->sc->comparison = SlepcCompareSmallestReal;
146: eps->sc->comparisonctx = NULL;
147: break;
148: case EPS_LARGEST_IMAGINARY:
149: eps->sc->comparison = SlepcCompareLargestImaginary;
150: eps->sc->comparisonctx = NULL;
151: break;
152: case EPS_SMALLEST_IMAGINARY:
153: eps->sc->comparison = SlepcCompareSmallestImaginary;
154: eps->sc->comparisonctx = NULL;
155: break;
156: case EPS_TARGET_MAGNITUDE:
157: eps->sc->comparison = SlepcCompareTargetMagnitude;
158: eps->sc->comparisonctx = &eps->target;
159: break;
160: case EPS_TARGET_REAL:
161: eps->sc->comparison = SlepcCompareTargetReal;
162: eps->sc->comparisonctx = &eps->target;
163: break;
164: case EPS_TARGET_IMAGINARY:
165: eps->sc->comparison = SlepcCompareTargetImaginary;
166: eps->sc->comparisonctx = &eps->target;
167: break;
168: case EPS_ALL:
169: eps->sc->comparison = SlepcCompareSmallestReal;
170: eps->sc->comparisonctx = NULL;
171: break;
172: case EPS_WHICH_USER:
173: break;
174: }
175: eps->sc->map = NULL;
176: eps->sc->mapobj = NULL;
178: /* fill sorting criterion for DS */
179: DSGetSlepcSC(eps->ds,&sc);
180: RGIsTrivial(eps->rg,&istrivial);
181: if (eps->which==EPS_ALL) {
182: sc->rg = NULL;
183: sc->comparison = SlepcCompareLargestMagnitude;
184: sc->comparisonctx = NULL;
185: sc->map = NULL;
186: sc->mapobj = NULL;
187: } else {
188: sc->rg = istrivial? NULL: eps->rg;
189: sc->comparison = eps->sc->comparison;
190: sc->comparisonctx = eps->sc->comparisonctx;
191: sc->map = SlepcMap_ST;
192: sc->mapobj = (PetscObject)eps->st;
193: }
195: /* Build balancing matrix if required */
196: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
197: if (!eps->D) {
198: BVCreateVec(eps->V,&eps->D);
199: PetscLogObjectParent((PetscObject)eps,(PetscObject)eps->D);
200: } else {
201: VecSet(eps->D,1.0);
202: }
203: EPSBuildBalance_Krylov(eps);
204: STSetBalanceMatrix(eps->st,eps->D);
205: }
207: /* Setup ST */
208: STSetUp(eps->st);
210: #if defined(PETSC_USE_COMPLEX)
211: STGetShift(eps->st,&sigma);
212: if (eps->ishermitian && PetscImaginaryPart(sigma) != 0.0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
213: #endif
214: PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg);
215: if (flg && eps->problem_type == EPS_PGNHEP) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
217: /* process deflation and initial vectors */
218: if (eps->nds<0) {
219: k = -eps->nds;
220: BVInsertConstraints(eps->V,&k,eps->defl);
221: SlepcBasisDestroy_Private(&eps->nds,&eps->defl);
222: eps->nds = k;
223: STCheckNullSpace(eps->st,eps->V);
224: }
225: if (eps->nini<0) {
226: k = -eps->nini;
227: if (k>eps->ncv) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The number of initial vectors is larger than ncv");
228: BVInsertVecs(eps->V,0,&k,eps->IS,PETSC_TRUE);
229: SlepcBasisDestroy_Private(&eps->nini,&eps->IS);
230: eps->nini = k;
231: }
233: PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
234: eps->state = EPS_STATE_SETUP;
235: return(0);
236: }
240: /*@
241: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
243: Collective on EPS and Mat
245: Input Parameters:
246: + eps - the eigenproblem solver context
247: . A - the matrix associated with the eigensystem
248: - B - the second matrix in the case of generalized eigenproblems
250: Notes:
251: To specify a standard eigenproblem, use NULL for parameter B.
253: It must be called before EPSSetUp(). If it is called again after EPSSetUp() then
254: the EPS object is reset.
256: Level: beginner
258: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetOperators()
259: @*/
260: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)261: {
263: PetscInt m,n,m0,nmat;
264: Mat mat[2];
273: /* Check for square matrices */
274: MatGetSize(A,&m,&n);
275: if (m!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix");
276: if (B) {
277: MatGetSize(B,&m0,&n);
278: if (m0!=n) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix");
279: if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match");
280: }
281: if (eps->state) { EPSReset(eps); }
282: eps->nrma = 0.0;
283: eps->nrmb = 0.0;
284: if (!eps->st) { EPSGetST(eps,&eps->st); }
285: mat[0] = A;
286: if (B) {
287: mat[1] = B;
288: nmat = 2;
289: } else nmat = 1;
290: STSetOperators(eps->st,nmat,mat);
291: return(0);
292: }
296: /*@
297: EPSGetOperators - Gets the matrices associated with the eigensystem.
299: Collective on EPS and Mat
301: Input Parameter:
302: . eps - the EPS context
304: Output Parameters:
305: + A - the matrix associated with the eigensystem
306: - B - the second matrix in the case of generalized eigenproblems
308: Level: intermediate
310: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
311: @*/
312: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)313: {
315: ST st;
316: PetscInt k;
320: EPSGetST(eps,&st);
321: if (A) { STGetOperators(st,0,A); }
322: if (B) {
323: STGetNumMatrices(st,&k);
324: if (k==1) B = NULL;
325: else {
326: STGetOperators(st,1,B);
327: }
328: }
329: return(0);
330: }
334: /*@
335: EPSSetDeflationSpace - Specify a basis of vectors that constitute the deflation
336: space.
338: Collective on EPS and Vec
340: Input Parameter:
341: + eps - the eigenproblem solver context
342: . n - number of vectors
343: - v - set of basis vectors of the deflation space
345: Notes:
346: When a deflation space is given, the eigensolver seeks the eigensolution
347: in the restriction of the problem to the orthogonal complement of this
348: space. This can be used for instance in the case that an invariant
349: subspace is known beforehand (such as the nullspace of the matrix).
351: These vectors do not persist from one EPSSolve() call to the other, so the
352: deflation space should be set every time.
354: The vectors do not need to be mutually orthonormal, since they are explicitly
355: orthonormalized internally.
357: Level: intermediate
359: .seealso: EPSSetInitialSpace()
360: @*/
361: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec *v)362: {
368: if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n out of range");
369: SlepcBasisReference_Private(n,v,&eps->nds,&eps->defl);
370: if (n>0) eps->state = EPS_STATE_INITIAL;
371: return(0);
372: }
376: /*@
377: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
378: space, that is, the subspace from which the solver starts to iterate.
380: Collective on EPS and Vec
382: Input Parameter:
383: + eps - the eigenproblem solver context
384: . n - number of vectors
385: - is - set of basis vectors of the initial space
387: Notes:
388: Some solvers start to iterate on a single vector (initial vector). In that case,
389: the other vectors are ignored.
391: These vectors do not persist from one EPSSolve() call to the other, so the
392: initial space should be set every time.
394: The vectors do not need to be mutually orthonormal, since they are explicitly
395: orthonormalized internally.
397: Common usage of this function is when the user can provide a rough approximation
398: of the wanted eigenspace. Then, convergence may be faster.
400: Level: intermediate
402: .seealso: EPSSetDeflationSpace()
403: @*/
404: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec *is)405: {
411: if (n<0) SETERRQ(PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
412: SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS);
413: if (n>0) eps->state = EPS_STATE_INITIAL;
414: return(0);
415: }
419: /*
420: EPSSetDimensions_Default - Set reasonable values for ncv, mpd if not set
421: by the user. This is called at setup.
422: */
423: PetscErrorCode EPSSetDimensions_Default(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)424: {
426: PetscBool krylov;
429: if (*ncv) { /* ncv set */
430: PetscObjectTypeCompareAny((PetscObject)eps,&krylov,EPSKRYLOVSCHUR,EPSARNOLDI,EPSLANCZOS,"");
431: if (krylov) {
432: if (*ncv<nev+1 && !(*ncv==nev && *ncv==eps->n)) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev+1");
433: } else {
434: if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)eps),1,"The value of ncv must be at least nev");
435: }
436: } else if (*mpd) { /* mpd set */
437: *ncv = PetscMin(eps->n,nev+(*mpd));
438: } else { /* neither set: defaults depend on nev being small or large */
439: if (nev<500) *ncv = PetscMin(eps->n,PetscMax(2*nev,nev+15));
440: else {
441: *mpd = 500;
442: *ncv = PetscMin(eps->n,nev+(*mpd));
443: }
444: }
445: if (!*mpd) *mpd = *ncv;
446: return(0);
447: }
451: /*@
452: EPSAllocateSolution - Allocate memory storage for common variables such
453: as eigenvalues and eigenvectors.
455: Collective on EPS457: Input Parameters:
458: + eps - eigensolver context
459: - extra - number of additional positions, used for methods that require a
460: working basis slightly larger than ncv
462: Developers Note:
463: This is PETSC_EXTERN because it may be required by user plugin EPS464: implementations.
466: Level: developer
467: @*/
468: PetscErrorCode EPSAllocateSolution(EPS eps,PetscInt extra)469: {
471: PetscInt oldsize,newc,requested;
472: PetscLogDouble cnt;
473: Vec t;
476: requested = eps->ncv + extra;
478: /* oldsize is zero if this is the first time setup is called */
479: BVGetSizes(eps->V,NULL,NULL,&oldsize);
480: newc = PetscMax(0,requested-oldsize);
482: /* allocate space for eigenvalues and friends */
483: if (requested != oldsize || !eps->eigr) {
484: if (oldsize) {
485: PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm);
486: }
487: PetscMalloc4(requested,&eps->eigr,requested,&eps->eigi,requested,&eps->errest,requested,&eps->perm);
488: cnt = 2*newc*sizeof(PetscScalar) + 2*newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
489: PetscLogObjectMemory((PetscObject)eps,cnt);
490: }
492: /* workspace for the case of arbitrary selection */
493: if (eps->arbitrary) {
494: if (eps->rr) {
495: PetscFree2(eps->rr,eps->ri);
496: }
497: PetscMalloc2(requested,&eps->rr,requested,&eps->ri);
498: PetscLogObjectMemory((PetscObject)eps,2*newc*sizeof(PetscScalar));
499: }
501: /* allocate V */
502: if (!eps->V) { EPSGetBV(eps,&eps->V); }
503: if (!oldsize) {
504: if (!((PetscObject)(eps->V))->type_name) {
505: BVSetType(eps->V,BVSVEC);
506: }
507: STMatCreateVecs(eps->st,&t,NULL);
508: BVSetSizesFromVec(eps->V,t,requested);
509: VecDestroy(&t);
510: } else {
511: BVResize(eps->V,requested,PETSC_FALSE);
512: }
513: return(0);
514: }