Actual source code: epssetup.c
slepc-3.20.1 2023-11-27
1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-, 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>
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: {
24: PetscFunctionBegin;
25: PetscTryTypeMethod(eps,setdefaultst);
26: if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSHIFT));
27: PetscFunctionReturn(PETSC_SUCCESS);
28: }
30: /*
31: This is done by preconditioned eigensolvers that use the PC only.
32: It sets STPRECOND with KSPPREONLY.
33: */
34: PetscErrorCode EPSSetDefaultST_Precond(EPS eps)
35: {
36: KSP ksp;
38: PetscFunctionBegin;
39: if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STPRECOND));
40: PetscCall(STGetKSP(eps->st,&ksp));
41: if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp,KSPPREONLY));
42: PetscFunctionReturn(PETSC_SUCCESS);
43: }
45: /*
46: This is done by preconditioned eigensolvers that can also use the KSP.
47: It sets STPRECOND with the default KSP (GMRES) and maxit=5.
48: */
49: PetscErrorCode EPSSetDefaultST_GMRES(EPS eps)
50: {
51: KSP ksp;
53: PetscFunctionBegin;
54: if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STPRECOND));
55: PetscCall(STPrecondSetKSPHasMat(eps->st,PETSC_TRUE));
56: PetscCall(STGetKSP(eps->st,&ksp));
57: if (!((PetscObject)ksp)->type_name) {
58: PetscCall(KSPSetType(ksp,KSPGMRES));
59: PetscCall(KSPSetTolerances(ksp,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT,5));
60: }
61: PetscFunctionReturn(PETSC_SUCCESS);
62: }
64: #if defined(SLEPC_HAVE_SCALAPACK) || defined(SLEPC_HAVE_ELPA) || defined(SLEPC_HAVE_ELEMENTAL) || defined(SLEPC_HAVE_EVSL)
65: /*
66: This is for direct eigensolvers that work with A and B directly, so
67: no need to factorize B.
68: */
69: PetscErrorCode EPSSetDefaultST_NoFactor(EPS eps)
70: {
71: KSP ksp;
72: PC pc;
74: PetscFunctionBegin;
75: if (!((PetscObject)eps->st)->type_name) PetscCall(STSetType(eps->st,STSHIFT));
76: PetscCall(STGetKSP(eps->st,&ksp));
77: if (!((PetscObject)ksp)->type_name) PetscCall(KSPSetType(ksp,KSPPREONLY));
78: PetscCall(KSPGetPC(ksp,&pc));
79: if (!((PetscObject)pc)->type_name) PetscCall(PCSetType(pc,PCNONE));
80: PetscFunctionReturn(PETSC_SUCCESS);
81: }
82: #endif
84: /*
85: Check that the ST selected by the user is compatible with the EPS solver and options
86: */
87: static PetscErrorCode EPSCheckCompatibleST(EPS eps)
88: {
89: PetscBool precond,shift,sinvert,cayley,lyapii;
90: #if defined(PETSC_USE_COMPLEX)
91: PetscScalar sigma;
92: #endif
94: PetscFunctionBegin;
95: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STPRECOND,&precond));
96: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSHIFT,&shift));
97: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STSINVERT,&sinvert));
98: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&cayley));
100: /* preconditioned eigensolvers */
101: PetscCheck(eps->categ!=EPS_CATEGORY_PRECOND || precond,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"This solver requires ST=PRECOND");
102: PetscCheck(eps->categ==EPS_CATEGORY_PRECOND || !precond,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"STPRECOND is intended for preconditioned eigensolvers only");
104: /* harmonic extraction */
105: PetscCheck(precond || shift || !eps->extraction || eps->extraction==EPS_RITZ,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cannot use a spectral transformation combined with harmonic extraction");
107: /* real shifts in Hermitian problems */
108: #if defined(PETSC_USE_COMPLEX)
109: PetscCall(STGetShift(eps->st,&sigma));
110: PetscCheck(!eps->ishermitian || PetscImaginaryPart(sigma)==0.0,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Hermitian problems are not compatible with complex shifts");
111: #endif
113: /* Cayley with PGNHEP */
114: PetscCheck(!cayley || eps->problem_type!=EPS_PGNHEP,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
116: /* make sure that the user does not specify smallest magnitude with shift-and-invert */
117: if ((cayley || sinvert) && (eps->categ==EPS_CATEGORY_KRYLOV || eps->categ==EPS_CATEGORY_OTHER)) {
118: PetscCall(PetscObjectTypeCompare((PetscObject)eps,EPSLYAPII,&lyapii));
119: PetscCheck(lyapii || eps->which==EPS_TARGET_MAGNITUDE || eps->which==EPS_TARGET_REAL || eps->which==EPS_TARGET_IMAGINARY || eps->which==EPS_ALL || eps->which==EPS_WHICH_USER,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"Shift-and-invert requires a target 'which' (see EPSSetWhichEigenpairs), for instance -st_type sinvert -eps_target 0 -eps_target_magnitude");
120: }
121: PetscFunctionReturn(PETSC_SUCCESS);
122: }
124: /*
125: MatEstimateSpectralRange_EPS: estimate the spectral range [left,right] of a
126: symmetric/Hermitian matrix A using an auxiliary EPS object
127: */
128: PetscErrorCode MatEstimateSpectralRange_EPS(Mat A,PetscReal *left,PetscReal *right)
129: {
130: PetscInt nconv;
131: PetscScalar eig0;
132: PetscReal tol=1e-3,errest=tol;
133: EPS eps;
135: PetscFunctionBegin;
136: *left = 0.0; *right = 0.0;
137: PetscCall(EPSCreate(PetscObjectComm((PetscObject)A),&eps));
138: PetscCall(EPSSetOptionsPrefix(eps,"eps_filter_"));
139: PetscCall(EPSSetOperators(eps,A,NULL));
140: PetscCall(EPSSetProblemType(eps,EPS_HEP));
141: PetscCall(EPSSetTolerances(eps,tol,50));
142: PetscCall(EPSSetConvergenceTest(eps,EPS_CONV_ABS));
143: PetscCall(EPSSetWhichEigenpairs(eps,EPS_SMALLEST_REAL));
144: PetscCall(EPSSolve(eps));
145: PetscCall(EPSGetConverged(eps,&nconv));
146: if (nconv>0) {
147: PetscCall(EPSGetEigenvalue(eps,0,&eig0,NULL));
148: PetscCall(EPSGetErrorEstimate(eps,0,&errest));
149: } else eig0 = eps->eigr[0];
150: *left = PetscRealPart(eig0)-errest;
151: PetscCall(EPSSetWhichEigenpairs(eps,EPS_LARGEST_REAL));
152: PetscCall(EPSSolve(eps));
153: PetscCall(EPSGetConverged(eps,&nconv));
154: if (nconv>0) {
155: PetscCall(EPSGetEigenvalue(eps,0,&eig0,NULL));
156: PetscCall(EPSGetErrorEstimate(eps,0,&errest));
157: } else eig0 = eps->eigr[0];
158: *right = PetscRealPart(eig0)+errest;
159: PetscCall(EPSDestroy(&eps));
160: PetscFunctionReturn(PETSC_SUCCESS);
161: }
163: /*
164: EPSSetUpSort_Basic: configure the EPS sorting criterion according to 'which'
165: */
166: PetscErrorCode EPSSetUpSort_Basic(EPS eps)
167: {
168: PetscFunctionBegin;
169: switch (eps->which) {
170: case EPS_LARGEST_MAGNITUDE:
171: eps->sc->comparison = SlepcCompareLargestMagnitude;
172: eps->sc->comparisonctx = NULL;
173: break;
174: case EPS_SMALLEST_MAGNITUDE:
175: eps->sc->comparison = SlepcCompareSmallestMagnitude;
176: eps->sc->comparisonctx = NULL;
177: break;
178: case EPS_LARGEST_REAL:
179: eps->sc->comparison = SlepcCompareLargestReal;
180: eps->sc->comparisonctx = NULL;
181: break;
182: case EPS_SMALLEST_REAL:
183: eps->sc->comparison = SlepcCompareSmallestReal;
184: eps->sc->comparisonctx = NULL;
185: break;
186: case EPS_LARGEST_IMAGINARY:
187: eps->sc->comparison = SlepcCompareLargestImaginary;
188: eps->sc->comparisonctx = NULL;
189: break;
190: case EPS_SMALLEST_IMAGINARY:
191: eps->sc->comparison = SlepcCompareSmallestImaginary;
192: eps->sc->comparisonctx = NULL;
193: break;
194: case EPS_TARGET_MAGNITUDE:
195: eps->sc->comparison = SlepcCompareTargetMagnitude;
196: eps->sc->comparisonctx = &eps->target;
197: break;
198: case EPS_TARGET_REAL:
199: eps->sc->comparison = SlepcCompareTargetReal;
200: eps->sc->comparisonctx = &eps->target;
201: break;
202: case EPS_TARGET_IMAGINARY:
203: #if defined(PETSC_USE_COMPLEX)
204: eps->sc->comparison = SlepcCompareTargetImaginary;
205: eps->sc->comparisonctx = &eps->target;
206: #endif
207: break;
208: case EPS_ALL:
209: eps->sc->comparison = SlepcCompareSmallestReal;
210: eps->sc->comparisonctx = NULL;
211: break;
212: case EPS_WHICH_USER:
213: break;
214: }
215: eps->sc->map = NULL;
216: eps->sc->mapobj = NULL;
217: PetscFunctionReturn(PETSC_SUCCESS);
218: }
220: /*
221: EPSSetUpSort_Default: configure both EPS and DS sorting criterion
222: */
223: PetscErrorCode EPSSetUpSort_Default(EPS eps)
224: {
225: SlepcSC sc;
226: PetscBool istrivial;
228: PetscFunctionBegin;
229: /* fill sorting criterion context */
230: PetscCall(EPSSetUpSort_Basic(eps));
231: /* fill sorting criterion for DS */
232: PetscCall(DSGetSlepcSC(eps->ds,&sc));
233: PetscCall(RGIsTrivial(eps->rg,&istrivial));
234: sc->rg = istrivial? NULL: eps->rg;
235: sc->comparison = eps->sc->comparison;
236: sc->comparisonctx = eps->sc->comparisonctx;
237: sc->map = SlepcMap_ST;
238: sc->mapobj = (PetscObject)eps->st;
239: PetscFunctionReturn(PETSC_SUCCESS);
240: }
242: /*@
243: EPSSetDSType - Sets the type of the internal DS object based on the current
244: settings of the eigenvalue solver.
246: Collective
248: Input Parameter:
249: . eps - eigenproblem solver context
251: Note:
252: This function need not be called explicitly, since it will be called at
253: both EPSSetFromOptions() and EPSSetUp().
255: Level: developer
257: .seealso: EPSSetFromOptions(), EPSSetUp()
258: @*/
259: PetscErrorCode EPSSetDSType(EPS eps)
260: {
261: PetscFunctionBegin;
263: PetscTryTypeMethod(eps,setdstype);
264: PetscFunctionReturn(PETSC_SUCCESS);
265: }
267: /*@
268: EPSSetUp - Sets up all the internal data structures necessary for the
269: execution of the eigensolver. Then calls STSetUp() for any set-up
270: operations associated to the ST object.
272: Collective
274: Input Parameter:
275: . eps - eigenproblem solver context
277: Notes:
278: This function need not be called explicitly in most cases, since EPSSolve()
279: calls it. It can be useful when one wants to measure the set-up time
280: separately from the solve time.
282: Level: developer
284: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp(), EPSSetInitialSpace()
285: @*/
286: PetscErrorCode EPSSetUp(EPS eps)
287: {
288: Mat A,B;
289: PetscInt k,nmat;
290: PetscBool flg;
292: PetscFunctionBegin;
294: if (eps->state) PetscFunctionReturn(PETSC_SUCCESS);
295: PetscCall(PetscLogEventBegin(EPS_SetUp,eps,0,0,0));
297: /* reset the convergence flag from the previous solves */
298: eps->reason = EPS_CONVERGED_ITERATING;
300: /* Set default solver type (EPSSetFromOptions was not called) */
301: if (!((PetscObject)eps)->type_name) PetscCall(EPSSetType(eps,EPSKRYLOVSCHUR));
302: if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
303: PetscCall(EPSSetDefaultST(eps));
305: PetscCall(STSetTransform(eps->st,PETSC_TRUE));
306: if (eps->useds && !eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
307: if (eps->useds) PetscCall(EPSSetDSType(eps));
308: if (eps->twosided) {
309: PetscCheck(!eps->ishermitian || (eps->isgeneralized && !eps->ispositive),PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Two-sided methods are not intended for %s problems",SLEPC_STRING_HERMITIAN);
310: }
311: if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
312: if (!((PetscObject)eps->rg)->type_name) PetscCall(RGSetType(eps->rg,RGINTERVAL));
314: /* Set problem dimensions */
315: PetscCall(STGetNumMatrices(eps->st,&nmat));
316: PetscCheck(nmat,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
317: PetscCall(STMatGetSize(eps->st,&eps->n,NULL));
318: PetscCall(STMatGetLocalSize(eps->st,&eps->nloc,NULL));
320: /* Set default problem type */
321: if (!eps->problem_type) {
322: if (nmat==1) PetscCall(EPSSetProblemType(eps,EPS_NHEP));
323: else PetscCall(EPSSetProblemType(eps,EPS_GNHEP));
324: } else if (nmat==1 && eps->isgeneralized) {
325: PetscCall(PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n"));
326: eps->isgeneralized = PETSC_FALSE;
327: eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
328: } else PetscCheck(nmat==1 || eps->isgeneralized,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Inconsistent EPS state: the problem type does not match the number of matrices");
330: if (eps->nev > eps->n) eps->nev = eps->n;
331: if (eps->ncv > eps->n) eps->ncv = eps->n;
333: /* check some combinations of eps->which */
334: PetscCheck(!eps->ishermitian || (eps->isgeneralized && !eps->ispositive) || (eps->which!=EPS_LARGEST_IMAGINARY && eps->which!=EPS_SMALLEST_IMAGINARY && eps->which!=EPS_TARGET_IMAGINARY),PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Sorting the eigenvalues along the imaginary axis is not allowed when all eigenvalues are real");
336: /* initialization of matrix norms */
337: if (eps->conv==EPS_CONV_NORM) {
338: if (!eps->nrma) {
339: PetscCall(STGetMatrix(eps->st,0,&A));
340: PetscCall(MatNorm(A,NORM_INFINITY,&eps->nrma));
341: }
342: if (nmat>1 && !eps->nrmb) {
343: PetscCall(STGetMatrix(eps->st,1,&B));
344: PetscCall(MatNorm(B,NORM_INFINITY,&eps->nrmb));
345: }
346: }
348: /* call specific solver setup */
349: PetscUseTypeMethod(eps,setup);
351: /* if purification is set, check that it really makes sense */
352: if (eps->purify) {
353: if (eps->categ==EPS_CATEGORY_PRECOND || eps->categ==EPS_CATEGORY_CONTOUR) eps->purify = PETSC_FALSE;
354: else {
355: if (!eps->isgeneralized) eps->purify = PETSC_FALSE;
356: else if (!eps->ishermitian && !eps->ispositive) eps->purify = PETSC_FALSE;
357: else {
358: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg));
359: if (flg) eps->purify = PETSC_FALSE;
360: }
361: }
362: }
364: /* set tolerance if not yet set */
365: if (eps->tol==(PetscReal)PETSC_DEFAULT) eps->tol = SLEPC_DEFAULT_TOL;
367: /* set up sorting criterion */
368: PetscTryTypeMethod(eps,setupsort);
370: /* Build balancing matrix if required */
371: if (eps->balance!=EPS_BALANCE_USER) {
372: PetscCall(STSetBalanceMatrix(eps->st,NULL));
373: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
374: if (!eps->D) PetscCall(BVCreateVec(eps->V,&eps->D));
375: PetscCall(EPSBuildBalance_Krylov(eps));
376: PetscCall(STSetBalanceMatrix(eps->st,eps->D));
377: }
378: }
380: /* Setup ST */
381: PetscCall(STSetUp(eps->st));
382: PetscCall(EPSCheckCompatibleST(eps));
384: /* process deflation and initial vectors */
385: if (eps->nds<0) {
386: k = -eps->nds;
387: PetscCall(BVInsertConstraints(eps->V,&k,eps->defl));
388: PetscCall(SlepcBasisDestroy_Private(&eps->nds,&eps->defl));
389: eps->nds = k;
390: PetscCall(STCheckNullSpace(eps->st,eps->V));
391: }
392: if (eps->nini<0) {
393: k = -eps->nini;
394: PetscCheck(k<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The number of initial vectors is larger than ncv");
395: PetscCall(BVInsertVecs(eps->V,0,&k,eps->IS,PETSC_TRUE));
396: PetscCall(SlepcBasisDestroy_Private(&eps->nini,&eps->IS));
397: eps->nini = k;
398: }
399: if (eps->twosided && eps->ninil<0) {
400: k = -eps->ninil;
401: PetscCheck(k<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The number of left initial vectors is larger than ncv");
402: PetscCall(BVInsertVecs(eps->W,0,&k,eps->ISL,PETSC_TRUE));
403: PetscCall(SlepcBasisDestroy_Private(&eps->ninil,&eps->ISL));
404: eps->ninil = k;
405: }
407: PetscCall(PetscLogEventEnd(EPS_SetUp,eps,0,0,0));
408: eps->state = EPS_STATE_SETUP;
409: PetscFunctionReturn(PETSC_SUCCESS);
410: }
412: /*@
413: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
415: Collective
417: Input Parameters:
418: + eps - the eigenproblem solver context
419: . A - the matrix associated with the eigensystem
420: - B - the second matrix in the case of generalized eigenproblems
422: Notes:
423: To specify a standard eigenproblem, use NULL for parameter B.
425: It must be called before EPSSetUp(). If it is called again after EPSSetUp() and
426: the matrix sizes have changed then the EPS object is reset.
428: Level: beginner
430: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetMatrix()
431: @*/
432: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
433: {
434: PetscInt m,n,m0,mloc,nloc,mloc0,nmat;
435: Mat mat[2];
437: PetscFunctionBegin;
441: PetscCheckSameComm(eps,1,A,2);
442: if (B) PetscCheckSameComm(eps,1,B,3);
444: /* Check matrix sizes */
445: PetscCall(MatGetSize(A,&m,&n));
446: PetscCall(MatGetLocalSize(A,&mloc,&nloc));
447: PetscCheck(m==n,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",m,n);
448: PetscCheck(mloc==nloc,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A does not have equal row and column sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc,nloc);
449: if (B) {
450: PetscCall(MatGetSize(B,&m0,&n));
451: PetscCall(MatGetLocalSize(B,&mloc0,&nloc));
452: PetscCheck(m0==n,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",m0,n);
453: PetscCheck(mloc0==nloc,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B does not have equal row and column local sizes (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc0,nloc);
454: PetscCheck(m==m0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match (%" PetscInt_FMT ", %" PetscInt_FMT ")",m,m0);
455: PetscCheck(mloc==mloc0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Local dimensions of A and B do not match (%" PetscInt_FMT ", %" PetscInt_FMT ")",mloc,mloc0);
456: }
457: if (eps->state && (n!=eps->n || nloc!=eps->nloc)) PetscCall(EPSReset(eps));
458: eps->nrma = 0.0;
459: eps->nrmb = 0.0;
460: if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
461: mat[0] = A;
462: if (B) {
463: mat[1] = B;
464: nmat = 2;
465: } else nmat = 1;
466: PetscCall(STSetMatrices(eps->st,nmat,mat));
467: eps->state = EPS_STATE_INITIAL;
468: PetscFunctionReturn(PETSC_SUCCESS);
469: }
471: /*@
472: EPSGetOperators - Gets the matrices associated with the eigensystem.
474: Collective
476: Input Parameter:
477: . eps - the EPS context
479: Output Parameters:
480: + A - the matrix associated with the eigensystem
481: - B - the second matrix in the case of generalized eigenproblems
483: Note:
484: Does not increase the reference count of the matrices, so you should not destroy them.
486: Level: intermediate
488: .seealso: EPSSolve(), EPSGetST(), STGetMatrix(), STSetMatrices()
489: @*/
490: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
491: {
492: ST st;
493: PetscInt k;
495: PetscFunctionBegin;
497: PetscCall(EPSGetST(eps,&st));
498: PetscCall(STGetNumMatrices(st,&k));
499: if (A) {
500: if (k<1) *A = NULL;
501: else PetscCall(STGetMatrix(st,0,A));
502: }
503: if (B) {
504: if (k<2) *B = NULL;
505: else PetscCall(STGetMatrix(st,1,B));
506: }
507: PetscFunctionReturn(PETSC_SUCCESS);
508: }
510: /*@C
511: EPSSetDeflationSpace - Specify a basis of vectors that constitute the deflation
512: space.
514: Collective
516: Input Parameters:
517: + eps - the eigenproblem solver context
518: . n - number of vectors
519: - v - set of basis vectors of the deflation space
521: Notes:
522: When a deflation space is given, the eigensolver seeks the eigensolution
523: in the restriction of the problem to the orthogonal complement of this
524: space. This can be used for instance in the case that an invariant
525: subspace is known beforehand (such as the nullspace of the matrix).
527: These vectors do not persist from one EPSSolve() call to the other, so the
528: deflation space should be set every time.
530: The vectors do not need to be mutually orthonormal, since they are explicitly
531: orthonormalized internally.
533: Level: intermediate
535: .seealso: EPSSetInitialSpace()
536: @*/
537: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec v[])
538: {
539: PetscFunctionBegin;
542: PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
543: if (n>0) {
544: PetscAssertPointer(v,3);
546: }
547: PetscCall(SlepcBasisReference_Private(n,v,&eps->nds,&eps->defl));
548: if (n>0) eps->state = EPS_STATE_INITIAL;
549: PetscFunctionReturn(PETSC_SUCCESS);
550: }
552: /*@C
553: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
554: space, that is, the subspace from which the solver starts to iterate.
556: Collective
558: Input Parameters:
559: + eps - the eigenproblem solver context
560: . n - number of vectors
561: - is - set of basis vectors of the initial space
563: Notes:
564: Some solvers start to iterate on a single vector (initial vector). In that case,
565: the other vectors are ignored.
567: These vectors do not persist from one EPSSolve() call to the other, so the
568: initial space should be set every time.
570: The vectors do not need to be mutually orthonormal, since they are explicitly
571: orthonormalized internally.
573: Common usage of this function is when the user can provide a rough approximation
574: of the wanted eigenspace. Then, convergence may be faster.
576: Level: intermediate
578: .seealso: EPSSetLeftInitialSpace(), EPSSetDeflationSpace()
579: @*/
580: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec is[])
581: {
582: PetscFunctionBegin;
585: PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
586: if (n>0) {
587: PetscAssertPointer(is,3);
589: }
590: PetscCall(SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS));
591: if (n>0) eps->state = EPS_STATE_INITIAL;
592: PetscFunctionReturn(PETSC_SUCCESS);
593: }
595: /*@C
596: EPSSetLeftInitialSpace - Specify a basis of vectors that constitute the left
597: initial space, used by two-sided solvers to start the left subspace.
599: Collective
601: Input Parameters:
602: + eps - the eigenproblem solver context
603: . n - number of vectors
604: - isl - set of basis vectors of the left initial space
606: Notes:
607: Left initial vectors are used to initiate the left search space in two-sided
608: eigensolvers. Users should pass here an approximation of the left eigenspace,
609: if available.
611: The same comments in EPSSetInitialSpace() are applicable here.
613: Level: intermediate
615: .seealso: EPSSetInitialSpace(), EPSSetTwoSided()
616: @*/
617: PetscErrorCode EPSSetLeftInitialSpace(EPS eps,PetscInt n,Vec isl[])
618: {
619: PetscFunctionBegin;
622: PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
623: if (n>0) {
624: PetscAssertPointer(isl,3);
626: }
627: PetscCall(SlepcBasisReference_Private(n,isl,&eps->ninil,&eps->ISL));
628: if (n>0) eps->state = EPS_STATE_INITIAL;
629: PetscFunctionReturn(PETSC_SUCCESS);
630: }
632: /*
633: EPSSetDimensions_Default - Set reasonable values for ncv, mpd if not set
634: by the user. This is called at setup.
635: */
636: PetscErrorCode EPSSetDimensions_Default(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
637: {
638: PetscBool krylov;
640: PetscFunctionBegin;
641: if (*ncv!=PETSC_DEFAULT) { /* ncv set */
642: PetscCall(PetscObjectTypeCompareAny((PetscObject)eps,&krylov,EPSKRYLOVSCHUR,EPSARNOLDI,EPSLANCZOS,""));
643: if (krylov) {
644: PetscCheck(*ncv>=nev+1 || (*ncv==nev && *ncv==eps->n),PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev+1");
645: } else {
646: PetscCheck(*ncv>=nev,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
647: }
648: } else if (*mpd!=PETSC_DEFAULT) { /* mpd set */
649: *ncv = PetscMin(eps->n,nev+(*mpd));
650: } else { /* neither set: defaults depend on nev being small or large */
651: if (nev<500) *ncv = PetscMin(eps->n,PetscMax(2*nev,nev+15));
652: else {
653: *mpd = 500;
654: *ncv = PetscMin(eps->n,nev+(*mpd));
655: }
656: }
657: if (*mpd==PETSC_DEFAULT) *mpd = *ncv;
658: PetscFunctionReturn(PETSC_SUCCESS);
659: }
661: /*@
662: EPSAllocateSolution - Allocate memory storage for common variables such
663: as eigenvalues and eigenvectors.
665: Collective
667: Input Parameters:
668: + eps - eigensolver context
669: - extra - number of additional positions, used for methods that require a
670: working basis slightly larger than ncv
672: Developer Notes:
673: This is SLEPC_EXTERN because it may be required by user plugin EPS
674: implementations.
676: Level: developer
678: .seealso: EPSSetUp()
679: @*/
680: PetscErrorCode EPSAllocateSolution(EPS eps,PetscInt extra)
681: {
682: PetscInt oldsize,requested;
683: PetscRandom rand;
684: Vec t;
686: PetscFunctionBegin;
687: requested = eps->ncv + extra;
689: /* oldsize is zero if this is the first time setup is called */
690: PetscCall(BVGetSizes(eps->V,NULL,NULL,&oldsize));
692: /* allocate space for eigenvalues and friends */
693: if (requested != oldsize || !eps->eigr) {
694: PetscCall(PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm));
695: PetscCall(PetscMalloc4(requested,&eps->eigr,requested,&eps->eigi,requested,&eps->errest,requested,&eps->perm));
696: }
698: /* workspace for the case of arbitrary selection */
699: if (eps->arbitrary) {
700: if (eps->rr) PetscCall(PetscFree2(eps->rr,eps->ri));
701: PetscCall(PetscMalloc2(requested,&eps->rr,requested,&eps->ri));
702: }
704: /* allocate V */
705: if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
706: if (!oldsize) {
707: if (!((PetscObject)(eps->V))->type_name) PetscCall(BVSetType(eps->V,BVMAT));
708: PetscCall(STMatCreateVecsEmpty(eps->st,&t,NULL));
709: PetscCall(BVSetSizesFromVec(eps->V,t,requested));
710: PetscCall(VecDestroy(&t));
711: } else PetscCall(BVResize(eps->V,requested,PETSC_FALSE));
713: /* allocate W */
714: if (eps->twosided) {
715: PetscCall(BVGetRandomContext(eps->V,&rand)); /* make sure the random context is available when duplicating */
716: PetscCall(BVDestroy(&eps->W));
717: PetscCall(BVDuplicate(eps->V,&eps->W));
718: }
719: PetscFunctionReturn(PETSC_SUCCESS);
720: }