Actual source code: epssetup.c
slepc-3.22.2 2024-12-02
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_CURRENT,PETSC_CURRENT,PETSC_CURRENT,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: PetscCall(STSetStructured(eps->st,PETSC_FALSE));
307: if (eps->useds && !eps->ds) PetscCall(EPSGetDS(eps,&eps->ds));
308: if (eps->useds) PetscCall(EPSSetDSType(eps));
309: if (eps->twosided) {
310: 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);
311: }
312: if (!eps->rg) PetscCall(EPSGetRG(eps,&eps->rg));
313: if (!((PetscObject)eps->rg)->type_name) PetscCall(RGSetType(eps->rg,RGINTERVAL));
315: /* Set problem dimensions */
316: PetscCall(STGetNumMatrices(eps->st,&nmat));
317: PetscCheck(nmat,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONGSTATE,"EPSSetOperators must be called first");
318: PetscCall(STMatGetSize(eps->st,&eps->n,NULL));
319: PetscCall(STMatGetLocalSize(eps->st,&eps->nloc,NULL));
321: /* Set default problem type */
322: if (!eps->problem_type) {
323: if (nmat==1) PetscCall(EPSSetProblemType(eps,EPS_NHEP));
324: else PetscCall(EPSSetProblemType(eps,EPS_GNHEP));
325: } else if (nmat==1 && eps->isgeneralized) {
326: PetscCall(PetscInfo(eps,"Eigenproblem set as generalized but no matrix B was provided; reverting to a standard eigenproblem\n"));
327: eps->isgeneralized = PETSC_FALSE;
328: eps->problem_type = eps->ishermitian? EPS_HEP: EPS_NHEP;
329: } 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");
331: if (eps->isstructured) {
332: /* make sure the user has set the appropriate matrix */
333: PetscCall(STGetMatrix(eps->st,0,&A));
334: if (eps->problem_type==EPS_BSE) PetscCall(SlepcCheckMatStruct(A,SLEPC_MAT_STRUCT_BSE,NULL));
335: }
337: /* safeguard for small problems */
338: if (eps->isstructured) {
339: if (2*eps->nev > eps->n) eps->nev = eps->n/2;
340: if (2*eps->ncv > eps->n) eps->ncv = eps->n/2;
341: } else {
342: if (eps->nev > eps->n) eps->nev = eps->n;
343: if (eps->ncv > eps->n) eps->ncv = eps->n;
344: }
346: /* check some combinations of eps->which */
347: 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");
349: /* initialization of matrix norms */
350: if (eps->conv==EPS_CONV_NORM) {
351: if (!eps->nrma) {
352: PetscCall(STGetMatrix(eps->st,0,&A));
353: PetscCall(MatNorm(A,NORM_INFINITY,&eps->nrma));
354: }
355: if (nmat>1 && !eps->nrmb) {
356: PetscCall(STGetMatrix(eps->st,1,&B));
357: PetscCall(MatNorm(B,NORM_INFINITY,&eps->nrmb));
358: }
359: }
361: /* call specific solver setup */
362: PetscUseTypeMethod(eps,setup);
364: /* if purification is set, check that it really makes sense */
365: if (eps->purify) {
366: if (eps->categ==EPS_CATEGORY_PRECOND || eps->categ==EPS_CATEGORY_CONTOUR) eps->purify = PETSC_FALSE;
367: else {
368: if (!eps->isgeneralized) eps->purify = PETSC_FALSE;
369: else if (!eps->ishermitian && !eps->ispositive) eps->purify = PETSC_FALSE;
370: else {
371: PetscCall(PetscObjectTypeCompare((PetscObject)eps->st,STCAYLEY,&flg));
372: if (flg) eps->purify = PETSC_FALSE;
373: }
374: }
375: }
377: /* set tolerance if not yet set */
378: if (eps->tol==(PetscReal)PETSC_DETERMINE) eps->tol = SLEPC_DEFAULT_TOL;
380: /* set up sorting criterion */
381: PetscTryTypeMethod(eps,setupsort);
383: /* Build balancing matrix if required */
384: if (eps->balance!=EPS_BALANCE_USER) {
385: PetscCall(STSetBalanceMatrix(eps->st,NULL));
386: if (!eps->ishermitian && (eps->balance==EPS_BALANCE_ONESIDE || eps->balance==EPS_BALANCE_TWOSIDE)) {
387: if (!eps->D) PetscCall(BVCreateVec(eps->V,&eps->D));
388: PetscCall(EPSBuildBalance_Krylov(eps));
389: PetscCall(STSetBalanceMatrix(eps->st,eps->D));
390: }
391: }
393: /* Setup ST */
394: PetscCall(STSetUp(eps->st));
395: PetscCall(EPSCheckCompatibleST(eps));
397: /* process deflation and initial vectors */
398: if (eps->nds<0) {
399: PetscCheck(!eps->isstructured,PetscObjectComm((PetscObject)eps),PETSC_ERR_SUP,"Deflation space is not supported in structured eigenproblems");
400: k = -eps->nds;
401: PetscCall(BVInsertConstraints(eps->V,&k,eps->defl));
402: PetscCall(SlepcBasisDestroy_Private(&eps->nds,&eps->defl));
403: eps->nds = k;
404: PetscCall(STCheckNullSpace(eps->st,eps->V));
405: }
406: if (eps->nini<0) {
407: k = -eps->nini;
408: PetscCheck(k<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The number of initial vectors is larger than ncv");
409: PetscCall(BVInsertVecs(eps->V,0,&k,eps->IS,PETSC_TRUE));
410: PetscCall(SlepcBasisDestroy_Private(&eps->nini,&eps->IS));
411: eps->nini = k;
412: }
413: if (eps->twosided && eps->ninil<0) {
414: k = -eps->ninil;
415: PetscCheck(k<=eps->ncv,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The number of left initial vectors is larger than ncv");
416: PetscCall(BVInsertVecs(eps->W,0,&k,eps->ISL,PETSC_TRUE));
417: PetscCall(SlepcBasisDestroy_Private(&eps->ninil,&eps->ISL));
418: eps->ninil = k;
419: }
421: PetscCall(PetscLogEventEnd(EPS_SetUp,eps,0,0,0));
422: eps->state = EPS_STATE_SETUP;
423: PetscFunctionReturn(PETSC_SUCCESS);
424: }
426: /*@
427: EPSSetOperators - Sets the matrices associated with the eigenvalue problem.
429: Collective
431: Input Parameters:
432: + eps - the eigenproblem solver context
433: . A - the matrix associated with the eigensystem
434: - B - the second matrix in the case of generalized eigenproblems
436: Notes:
437: To specify a standard eigenproblem, use NULL for parameter B.
439: It must be called before EPSSetUp(). If it is called again after EPSSetUp() and
440: the matrix sizes have changed then the EPS object is reset.
442: For structured eigenproblem types such as EPS_BSE (see EPSSetProblemType()), the
443: provided matrices must have been created with the corresponding helper function,
444: i.e., MatCreateBSE().
446: Level: beginner
448: .seealso: EPSSolve(), EPSSetUp(), EPSReset(), EPSGetST(), STGetMatrix(), EPSSetProblemType()
449: @*/
450: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
451: {
452: PetscInt m,n,m0,mloc,nloc,mloc0,nmat;
453: Mat mat[2];
455: PetscFunctionBegin;
459: PetscCheckSameComm(eps,1,A,2);
460: if (B) PetscCheckSameComm(eps,1,B,3);
462: /* Check matrix sizes */
463: PetscCall(MatGetSize(A,&m,&n));
464: PetscCall(MatGetLocalSize(A,&mloc,&nloc));
465: PetscCheck(m==n,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"A is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",m,n);
466: 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);
467: if (B) {
468: PetscCall(MatGetSize(B,&m0,&n));
469: PetscCall(MatGetLocalSize(B,&mloc0,&nloc));
470: PetscCheck(m0==n,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_WRONG,"B is a non-square matrix (%" PetscInt_FMT " rows, %" PetscInt_FMT " cols)",m0,n);
471: 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);
472: PetscCheck(m==m0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_INCOMP,"Dimensions of A and B do not match (%" PetscInt_FMT ", %" PetscInt_FMT ")",m,m0);
473: 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);
474: }
475: if (eps->state && (n!=eps->n || nloc!=eps->nloc)) PetscCall(EPSReset(eps));
476: eps->nrma = 0.0;
477: eps->nrmb = 0.0;
478: if (!eps->st) PetscCall(EPSGetST(eps,&eps->st));
479: mat[0] = A;
480: if (B) {
481: mat[1] = B;
482: nmat = 2;
483: } else nmat = 1;
484: PetscCall(STSetMatrices(eps->st,nmat,mat));
485: eps->state = EPS_STATE_INITIAL;
486: PetscFunctionReturn(PETSC_SUCCESS);
487: }
489: /*@
490: EPSGetOperators - Gets the matrices associated with the eigensystem.
492: Collective
494: Input Parameter:
495: . eps - the EPS context
497: Output Parameters:
498: + A - the matrix associated with the eigensystem
499: - B - the second matrix in the case of generalized eigenproblems
501: Note:
502: Does not increase the reference count of the matrices, so you should not destroy them.
504: Level: intermediate
506: .seealso: EPSSolve(), EPSGetST(), STGetMatrix(), STSetMatrices()
507: @*/
508: PetscErrorCode EPSGetOperators(EPS eps,Mat *A,Mat *B)
509: {
510: ST st;
511: PetscInt k;
513: PetscFunctionBegin;
515: PetscCall(EPSGetST(eps,&st));
516: PetscCall(STGetNumMatrices(st,&k));
517: if (A) {
518: if (k<1) *A = NULL;
519: else PetscCall(STGetMatrix(st,0,A));
520: }
521: if (B) {
522: if (k<2) *B = NULL;
523: else PetscCall(STGetMatrix(st,1,B));
524: }
525: PetscFunctionReturn(PETSC_SUCCESS);
526: }
528: /*@
529: EPSSetDeflationSpace - Specify a basis of vectors that constitute the deflation
530: space.
532: Collective
534: Input Parameters:
535: + eps - the eigenproblem solver context
536: . n - number of vectors
537: - v - set of basis vectors of the deflation space
539: Notes:
540: When a deflation space is given, the eigensolver seeks the eigensolution
541: in the restriction of the problem to the orthogonal complement of this
542: space. This can be used for instance in the case that an invariant
543: subspace is known beforehand (such as the nullspace of the matrix).
545: These vectors do not persist from one EPSSolve() call to the other, so the
546: deflation space should be set every time.
548: The vectors do not need to be mutually orthonormal, since they are explicitly
549: orthonormalized internally.
551: Level: intermediate
553: .seealso: EPSSetInitialSpace()
554: @*/
555: PetscErrorCode EPSSetDeflationSpace(EPS eps,PetscInt n,Vec v[])
556: {
557: PetscFunctionBegin;
560: PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
561: if (n>0) {
562: PetscAssertPointer(v,3);
564: }
565: PetscCall(SlepcBasisReference_Private(n,v,&eps->nds,&eps->defl));
566: if (n>0) eps->state = EPS_STATE_INITIAL;
567: PetscFunctionReturn(PETSC_SUCCESS);
568: }
570: /*@
571: EPSSetInitialSpace - Specify a basis of vectors that constitute the initial
572: space, that is, the subspace from which the solver starts to iterate.
574: Collective
576: Input Parameters:
577: + eps - the eigenproblem solver context
578: . n - number of vectors
579: - is - set of basis vectors of the initial space
581: Notes:
582: Some solvers start to iterate on a single vector (initial vector). In that case,
583: the other vectors are ignored.
585: These vectors do not persist from one EPSSolve() call to the other, so the
586: initial space should be set every time.
588: The vectors do not need to be mutually orthonormal, since they are explicitly
589: orthonormalized internally.
591: Common usage of this function is when the user can provide a rough approximation
592: of the wanted eigenspace. Then, convergence may be faster.
594: Level: intermediate
596: .seealso: EPSSetLeftInitialSpace(), EPSSetDeflationSpace()
597: @*/
598: PetscErrorCode EPSSetInitialSpace(EPS eps,PetscInt n,Vec is[])
599: {
600: PetscFunctionBegin;
603: PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
604: if (n>0) {
605: PetscAssertPointer(is,3);
607: }
608: PetscCall(SlepcBasisReference_Private(n,is,&eps->nini,&eps->IS));
609: if (n>0) eps->state = EPS_STATE_INITIAL;
610: PetscFunctionReturn(PETSC_SUCCESS);
611: }
613: /*@
614: EPSSetLeftInitialSpace - Specify a basis of vectors that constitute the left
615: initial space, used by two-sided solvers to start the left subspace.
617: Collective
619: Input Parameters:
620: + eps - the eigenproblem solver context
621: . n - number of vectors
622: - isl - set of basis vectors of the left initial space
624: Notes:
625: Left initial vectors are used to initiate the left search space in two-sided
626: eigensolvers. Users should pass here an approximation of the left eigenspace,
627: if available.
629: The same comments in EPSSetInitialSpace() are applicable here.
631: Level: intermediate
633: .seealso: EPSSetInitialSpace(), EPSSetTwoSided()
634: @*/
635: PetscErrorCode EPSSetLeftInitialSpace(EPS eps,PetscInt n,Vec isl[])
636: {
637: PetscFunctionBegin;
640: PetscCheck(n>=0,PetscObjectComm((PetscObject)eps),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
641: if (n>0) {
642: PetscAssertPointer(isl,3);
644: }
645: PetscCall(SlepcBasisReference_Private(n,isl,&eps->ninil,&eps->ISL));
646: if (n>0) eps->state = EPS_STATE_INITIAL;
647: PetscFunctionReturn(PETSC_SUCCESS);
648: }
650: /*
651: EPSSetDimensions_Default - Set reasonable values for ncv, mpd if not set
652: by the user. This is called at setup.
653: */
654: PetscErrorCode EPSSetDimensions_Default(EPS eps,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
655: {
656: PetscBool krylov;
657: PetscInt n = eps->isstructured? eps->n/2: eps->n;
659: PetscFunctionBegin;
660: if (*ncv!=PETSC_DETERMINE) { /* ncv set */
661: PetscCall(PetscObjectTypeCompareAny((PetscObject)eps,&krylov,EPSKRYLOVSCHUR,EPSARNOLDI,EPSLANCZOS,""));
662: if (krylov) {
663: PetscCheck(*ncv>=nev+1 || (*ncv==nev && *ncv==n),PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev+1");
664: } else {
665: PetscCheck(*ncv>=nev,PetscObjectComm((PetscObject)eps),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
666: }
667: } else if (*mpd!=PETSC_DETERMINE) { /* mpd set */
668: *ncv = PetscMin(n,nev+(*mpd));
669: } else { /* neither set: defaults depend on nev being small or large */
670: if (nev<500) *ncv = PetscMin(n,PetscMax(2*nev,nev+15));
671: else {
672: *mpd = 500;
673: *ncv = PetscMin(n,nev+(*mpd));
674: }
675: }
676: if (*mpd==PETSC_DETERMINE) *mpd = *ncv;
677: PetscFunctionReturn(PETSC_SUCCESS);
678: }
680: /*@
681: EPSAllocateSolution - Allocate memory storage for common variables such
682: as eigenvalues and eigenvectors.
684: Collective
686: Input Parameters:
687: + eps - eigensolver context
688: - extra - number of additional positions, used for methods that require a
689: working basis slightly larger than ncv
691: Developer Notes:
692: This is SLEPC_EXTERN because it may be required by user plugin EPS
693: implementations.
695: Level: developer
697: .seealso: EPSSetUp()
698: @*/
699: PetscErrorCode EPSAllocateSolution(EPS eps,PetscInt extra)
700: {
701: PetscInt oldsize,requested;
702: PetscRandom rand;
703: Vec t;
705: PetscFunctionBegin;
706: requested = eps->ncv + extra;
708: /* oldsize is zero if this is the first time setup is called */
709: PetscCall(BVGetSizes(eps->V,NULL,NULL,&oldsize));
711: /* allocate space for eigenvalues and friends */
712: if (requested != oldsize || !eps->eigr) {
713: PetscCall(PetscFree4(eps->eigr,eps->eigi,eps->errest,eps->perm));
714: PetscCall(PetscMalloc4(requested,&eps->eigr,requested,&eps->eigi,requested,&eps->errest,requested,&eps->perm));
715: }
717: /* workspace for the case of arbitrary selection */
718: if (eps->arbitrary) {
719: if (eps->rr) PetscCall(PetscFree2(eps->rr,eps->ri));
720: PetscCall(PetscMalloc2(requested,&eps->rr,requested,&eps->ri));
721: }
723: /* allocate V */
724: if (!eps->V) PetscCall(EPSGetBV(eps,&eps->V));
725: if (!oldsize) {
726: if (!((PetscObject)eps->V)->type_name) PetscCall(BVSetType(eps->V,BVMAT));
727: PetscCall(STMatCreateVecsEmpty(eps->st,&t,NULL));
728: PetscCall(BVSetSizesFromVec(eps->V,t,requested));
729: PetscCall(VecDestroy(&t));
730: } else PetscCall(BVResize(eps->V,requested,PETSC_FALSE));
732: /* allocate W */
733: if (eps->twosided) {
734: PetscCall(BVGetRandomContext(eps->V,&rand)); /* make sure the random context is available when duplicating */
735: PetscCall(BVDestroy(&eps->W));
736: PetscCall(BVDuplicate(eps->V,&eps->W));
737: }
738: PetscFunctionReturn(PETSC_SUCCESS);
739: }