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: PEP routines related to problem setup
12: */
14: #include <slepc/private/pepimpl.h> /*I "slepcpep.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 PEPSetFromOptions (before STSetFromOptions)
20: and also at PEPSetUp (in case PEPSetFromOptions was not called).
21: */
22: PetscErrorCode PEPSetDefaultST(PEP pep) 23: {
27: if (pep->ops->setdefaultst) { (*pep->ops->setdefaultst)(pep); }
28: if (!((PetscObject)pep->st)->type_name) {
29: STSetType(pep->st,STSHIFT);
30: }
31: return(0);
32: }
34: /*
35: This is used in Q-Arnoldi and STOAR to set the transform flag by
36: default, otherwise the user has to explicitly run with -st_transform
37: */
38: PetscErrorCode PEPSetDefaultST_Transform(PEP pep) 39: {
43: STSetTransform(pep->st,PETSC_TRUE);
44: return(0);
45: }
47: /*@
48: PEPSetUp - Sets up all the internal data structures necessary for the
49: execution of the PEP solver.
51: Collective on PEP 53: Input Parameter:
54: . pep - solver context
56: Notes:
57: This function need not be called explicitly in most cases, since PEPSolve()
58: calls it. It can be useful when one wants to measure the set-up time
59: separately from the solve time.
61: Level: developer
63: .seealso: PEPCreate(), PEPSolve(), PEPDestroy()
64: @*/
65: PetscErrorCode PEPSetUp(PEP pep) 66: {
68: SlepcSC sc;
69: PetscBool istrivial,flg;
70: PetscInt k;
71: KSP ksp;
72: PC pc;
73: PetscMPIInt size;
74: const MatSolverPackage stype;
78: if (pep->state) return(0);
79: PetscLogEventBegin(PEP_SetUp,pep,0,0,0);
81: /* reset the convergence flag from the previous solves */
82: pep->reason = PEP_CONVERGED_ITERATING;
84: /* set default solver type (PEPSetFromOptions was not called) */
85: if (!((PetscObject)pep)->type_name) {
86: PEPSetType(pep,PEPTOAR);
87: }
88: if (!pep->st) { PEPGetST(pep,&pep->st); }
89: PEPSetDefaultST(pep);
90: if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
91: if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
92: if (!((PetscObject)pep->rg)->type_name) {
93: RGSetType(pep->rg,RGINTERVAL);
94: }
96: /* check matrices, transfer them to ST */
97: if (!pep->A) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"PEPSetOperators must be called first");
98: STSetMatrices(pep->st,pep->nmat,pep->A);
100: /* set problem dimensions */
101: MatGetSize(pep->A[0],&pep->n,NULL);
102: MatGetLocalSize(pep->A[0],&pep->nloc,NULL);
104: /* set default problem type */
105: if (!pep->problem_type) {
106: PEPSetProblemType(pep,PEP_GENERAL);
107: }
109: /* check consistency of refinement options */
110: if (pep->refine) {
111: if (!pep->scheme) { /* set default scheme */
112: PEPRefineGetKSP(pep,&ksp);
113: KSPGetPC(ksp,&pc);
114: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
115: if (flg) {
116: PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
117: }
118: pep->scheme = flg? PEP_REFINE_SCHEME_MBE: PEP_REFINE_SCHEME_SCHUR;
119: }
120: if (pep->scheme==PEP_REFINE_SCHEME_MBE) {
121: PEPRefineGetKSP(pep,&ksp);
122: KSPGetPC(ksp,&pc);
123: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
124: if (flg) {
125: PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
126: }
127: if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"The MBE scheme for refinement requires a direct solver in KSP");
128: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
129: if (size>1) { /* currently selected PC is a factorization */
130: PCFactorGetMatSolverPackage(pc,&stype);
131: PetscStrcmp(stype,MATSOLVERPETSC,&flg);
132: if (flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"For Newton refinement, you chose to solve linear systems with a factorization, but in parallel runs you need to select an external package");
133: }
134: }
135: if (pep->scheme==PEP_REFINE_SCHEME_SCHUR) {
136: if (pep->npart>1) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"The Schur scheme for refinement does not support subcommunicators");
137: }
138: }
139: /* call specific solver setup */
140: (*pep->ops->setup)(pep);
142: /* set tolerance if not yet set */
143: if (pep->tol==PETSC_DEFAULT) pep->tol = SLEPC_DEFAULT_TOL;
144: if (pep->refine) {
145: if (pep->rtol==PETSC_DEFAULT) pep->rtol = PetscMax(pep->tol/1000,PETSC_MACHINE_EPSILON);
146: if (pep->rits==PETSC_DEFAULT) pep->rits = (pep->refine==PEP_REFINE_SIMPLE)? 10: 1;
147: }
149: /* set default extraction */
150: if (!pep->extract) {
151: pep->extract = (pep->basis==PEP_BASIS_MONOMIAL)? PEP_EXTRACT_NORM: PEP_EXTRACT_NONE;
152: }
154: /* fill sorting criterion context */
155: switch (pep->which) {
156: case PEP_LARGEST_MAGNITUDE:
157: pep->sc->comparison = SlepcCompareLargestMagnitude;
158: pep->sc->comparisonctx = NULL;
159: break;
160: case PEP_SMALLEST_MAGNITUDE:
161: pep->sc->comparison = SlepcCompareSmallestMagnitude;
162: pep->sc->comparisonctx = NULL;
163: break;
164: case PEP_LARGEST_REAL:
165: pep->sc->comparison = SlepcCompareLargestReal;
166: pep->sc->comparisonctx = NULL;
167: break;
168: case PEP_SMALLEST_REAL:
169: pep->sc->comparison = SlepcCompareSmallestReal;
170: pep->sc->comparisonctx = NULL;
171: break;
172: case PEP_LARGEST_IMAGINARY:
173: pep->sc->comparison = SlepcCompareLargestImaginary;
174: pep->sc->comparisonctx = NULL;
175: break;
176: case PEP_SMALLEST_IMAGINARY:
177: pep->sc->comparison = SlepcCompareSmallestImaginary;
178: pep->sc->comparisonctx = NULL;
179: break;
180: case PEP_TARGET_MAGNITUDE:
181: pep->sc->comparison = SlepcCompareTargetMagnitude;
182: pep->sc->comparisonctx = &pep->target;
183: break;
184: case PEP_TARGET_REAL:
185: pep->sc->comparison = SlepcCompareTargetReal;
186: pep->sc->comparisonctx = &pep->target;
187: break;
188: case PEP_TARGET_IMAGINARY:
189: #if defined(PETSC_USE_COMPLEX)
190: pep->sc->comparison = SlepcCompareTargetImaginary;
191: pep->sc->comparisonctx = &pep->target;
192: #endif
193: break;
194: case PEP_WHICH_USER:
195: break;
196: }
197: pep->sc->map = NULL;
198: pep->sc->mapobj = NULL;
200: /* fill sorting criterion for DS */
201: DSGetSlepcSC(pep->ds,&sc);
202: RGIsTrivial(pep->rg,&istrivial);
203: sc->rg = istrivial? NULL: pep->rg;
204: sc->comparison = pep->sc->comparison;
205: sc->comparisonctx = pep->sc->comparisonctx;
206: sc->map = SlepcMap_ST;
207: sc->mapobj = (PetscObject)pep->st;
209: /* setup ST */
210: STSetUp(pep->st);
212: /* compute matrix coefficients */
213: STGetTransform(pep->st,&flg);
214: if (!flg) {
215: if (pep->solvematcoeffs) { STMatSetUp(pep->st,1.0,pep->solvematcoeffs); }
216: } else {
217: if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Cannot use ST-transform with non-monomial basis in PEP");
218: }
220: /* compute scale factor if no set by user */
221: PEPComputeScaleFactor(pep);
223: /* build balancing matrix if required */
224: if (pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) {
225: if (!pep->Dl) {
226: BVCreateVec(pep->V,&pep->Dl);
227: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->Dl);
228: }
229: if (!pep->Dr) {
230: BVCreateVec(pep->V,&pep->Dr);
231: PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->Dr);
232: }
233: PEPBuildDiagonalScaling(pep);
234: }
236: /* process initial vectors */
237: if (pep->nini<0) {
238: k = -pep->nini;
239: if (k>pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),1,"The number of initial vectors is larger than ncv");
240: BVInsertVecs(pep->V,0,&k,pep->IS,PETSC_TRUE);
241: SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
242: pep->nini = k;
243: }
244: PetscLogEventEnd(PEP_SetUp,pep,0,0,0);
245: pep->state = PEP_STATE_SETUP;
246: return(0);
247: }
249: /*@
250: PEPSetOperators - Sets the coefficient matrices associated with the polynomial
251: eigenvalue problem.
253: Collective on PEP and Mat
255: Input Parameters:
256: + pep - the eigenproblem solver context
257: . nmat - number of matrices in array A
258: - A - the array of matrices associated with the eigenproblem
260: Notes:
261: The polynomial eigenproblem is defined as P(l)*x=0, where l is
262: the eigenvalue, x is the eigenvector, and P(l) is defined as
263: P(l) = A_0 + l*A_1 + ... + l^d*A_d, with d=nmat-1 (the degree of P).
264: For non-monomial bases, this expression is different.
266: Level: beginner
268: .seealso: PEPSolve(), PEPGetOperators(), PEPGetNumMatrices(), PEPSetBasis()
269: @*/
270: PetscErrorCode PEPSetOperators(PEP pep,PetscInt nmat,Mat A[])271: {
273: PetscInt i,n,m,m0=0;
278: if (nmat <= 0) SETERRQ1(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Non-positive value of nmat: %D",nmat);
279: if (nmat <= 2) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Cannot solve linear eigenproblems with PEP; use EPS instead");
283: MatGetSize(A[0],&m,&n);
284: if (pep->state && n!=pep->n) { PEPReset(pep); }
285: else if (pep->nmat) {
286: MatDestroyMatrices(pep->nmat,&pep->A);
287: PetscFree2(pep->pbc,pep->nrma);
288: PetscFree(pep->solvematcoeffs);
289: }
291: PetscMalloc1(nmat,&pep->A);
292: PetscCalloc2(3*nmat,&pep->pbc,nmat,&pep->nrma);
293: for (i=0;i<nmat;i++) pep->pbc[i] = 1.0; /* default to monomial basis */
294: PetscLogObjectMemory((PetscObject)pep,nmat*sizeof(Mat)+4*nmat*sizeof(PetscReal));
295: for (i=0;i<nmat;i++) {
298: MatGetSize(A[i],&m,&n);
299: if (m!=n) SETERRQ1(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"A[%D] is a non-square matrix",i);
300: if (!i) m0 = m;
301: if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_INCOMP,"Dimensions of matrices do not match with each other");
302: PetscObjectReference((PetscObject)A[i]);
303: pep->A[i] = A[i];
304: }
305: pep->nmat = nmat;
306: pep->state = PEP_STATE_INITIAL;
307: return(0);
308: }
310: /*@
311: PEPGetOperators - Gets the matrices associated with the polynomial eigensystem.
313: Not collective, though parallel Mats are returned if the PEP is parallel
315: Input Parameters:
316: + pep - the PEP context
317: - k - the index of the requested matrix (starting in 0)
319: Output Parameter:
320: . A - the requested matrix
322: Level: intermediate
324: .seealso: PEPSolve(), PEPSetOperators(), PEPGetNumMatrices()
325: @*/
326: PetscErrorCode PEPGetOperators(PEP pep,PetscInt k,Mat *A)327: {
331: if (k<0 || k>=pep->nmat) SETERRQ1(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",pep->nmat-1);
332: *A = pep->A[k];
333: return(0);
334: }
336: /*@
337: PEPGetNumMatrices - Returns the number of matrices stored in the PEP.
339: Not collective
341: Input Parameter:
342: . pep - the PEP context
344: Output Parameters:
345: . nmat - the number of matrices passed in PEPSetOperators()
347: Level: intermediate
349: .seealso: PEPSetOperators()
350: @*/
351: PetscErrorCode PEPGetNumMatrices(PEP pep,PetscInt *nmat)352: {
356: *nmat = pep->nmat;
357: return(0);
358: }
360: /*@
361: PEPSetInitialSpace - Specify a basis of vectors that constitute the initial
362: space, that is, the subspace from which the solver starts to iterate.
364: Collective on PEP and Vec
366: Input Parameter:
367: + pep - the polynomial eigensolver context
368: . n - number of vectors
369: - is - set of basis vectors of the initial space
371: Notes:
372: Some solvers start to iterate on a single vector (initial vector). In that case,
373: the other vectors are ignored.
375: These vectors do not persist from one PEPSolve() call to the other, so the
376: initial space should be set every time.
378: The vectors do not need to be mutually orthonormal, since they are explicitly
379: orthonormalized internally.
381: Common usage of this function is when the user can provide a rough approximation
382: of the wanted eigenspace. Then, convergence may be faster.
384: Level: intermediate
385: @*/
386: PetscErrorCode PEPSetInitialSpace(PEP pep,PetscInt n,Vec *is)387: {
393: if (n<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
394: SlepcBasisReference_Private(n,is,&pep->nini,&pep->IS);
395: if (n>0) pep->state = PEP_STATE_INITIAL;
396: return(0);
397: }
399: /*
400: PEPSetDimensions_Default - Set reasonable values for ncv, mpd if not set
401: by the user. This is called at setup.
402: */
403: PetscErrorCode PEPSetDimensions_Default(PEP pep,PetscInt nev,PetscInt *ncv,PetscInt *mpd)404: {
406: PetscBool krylov;
407: PetscInt dim;
410: PetscObjectTypeCompareAny((PetscObject)pep,&krylov,PEPTOAR,PEPQARNOLDI,"");
411: dim = krylov?(pep->nmat-1)*pep->n:pep->n;
412: if (*ncv) { /* ncv set */
413: if (krylov) {
414: if (*ncv<nev+1 && !(*ncv==nev && *ncv==dim)) SETERRQ(PetscObjectComm((PetscObject)pep),1,"The value of ncv must be at least nev+1");
415: } else {
416: if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)pep),1,"The value of ncv must be at least nev");
417: }
418: } else if (*mpd) { /* mpd set */
419: *ncv = PetscMin(dim,nev+(*mpd));
420: } else { /* neither set: defaults depend on nev being small or large */
421: if (nev<500) *ncv = PetscMin(dim,PetscMax(2*nev,nev+15));
422: else {
423: *mpd = 500;
424: *ncv = PetscMin(dim,nev+(*mpd));
425: }
426: }
427: if (!*mpd) *mpd = *ncv;
428: return(0);
429: }
431: /*@
432: PEPAllocateSolution - Allocate memory storage for common variables such
433: as eigenvalues and eigenvectors.
435: Collective on PEP437: Input Parameters:
438: + pep - eigensolver context
439: - extra - number of additional positions, used for methods that require a
440: working basis slightly larger than ncv
442: Developers Note:
443: This is PETSC_EXTERN because it may be required by user plugin PEP444: implementations.
446: Level: developer
447: @*/
448: PetscErrorCode PEPAllocateSolution(PEP pep,PetscInt extra)449: {
451: PetscInt oldsize,newc,requested,requestedbv;
452: PetscLogDouble cnt;
453: Vec t;
456: requested = (pep->lineariz? pep->ncv: pep->ncv*(pep->nmat-1)) + extra;
457: requestedbv = pep->ncv + extra;
459: /* oldsize is zero if this is the first time setup is called */
460: BVGetSizes(pep->V,NULL,NULL,&oldsize);
462: /* allocate space for eigenvalues and friends */
463: if (requested != oldsize || !pep->eigr) {
464: PetscFree4(pep->eigr,pep->eigi,pep->errest,pep->perm);
465: PetscMalloc4(requested,&pep->eigr,requested,&pep->eigi,requested,&pep->errest,requested,&pep->perm);
466: newc = PetscMax(0,requested-oldsize);
467: cnt = 2*newc*sizeof(PetscScalar) + newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
468: PetscLogObjectMemory((PetscObject)pep,cnt);
469: }
471: /* allocate V */
472: if (!pep->V) { PEPGetBV(pep,&pep->V); }
473: if (!oldsize) {
474: if (!((PetscObject)(pep->V))->type_name) {
475: BVSetType(pep->V,BVSVEC);
476: }
477: STMatCreateVecsEmpty(pep->st,&t,NULL);
478: BVSetSizesFromVec(pep->V,t,requestedbv);
479: VecDestroy(&t);
480: } else {
481: BVResize(pep->V,requestedbv,PETSC_FALSE);
482: }
483: return(0);
484: }