1: /*
2: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
3: SLEPc - Scalable Library for Eigenvalue Problem Computations
4: Copyright (c) 2002-2020, 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: NEP routines related to problem setup
12: */
14: #include <slepc/private/nepimpl.h> /*I "slepcnep.h" I*/
16: /*@
17: NEPSetUp - Sets up all the internal data structures necessary for the
18: execution of the NEP solver.
20: Collective on nep
22: Input Parameter:
23: . nep - solver context
25: Notes:
26: This function need not be called explicitly in most cases, since NEPSolve()
27: calls it. It can be useful when one wants to measure the set-up time
28: separately from the solve time.
30: Level: developer
32: .seealso: NEPCreate(), NEPSolve(), NEPDestroy()
33: @*/
34: PetscErrorCode NEPSetUp(NEP nep) 35: {
37: PetscInt k;
38: SlepcSC sc;
39: Mat T;
40: PetscBool flg;
41: KSP ksp;
42: PC pc;
43: PetscMPIInt size;
44: MatSolverType stype;
48: NEPCheckProblem(nep,1);
49: if (nep->state) return(0);
50: PetscLogEventBegin(NEP_SetUp,nep,0,0,0);
52: /* reset the convergence flag from the previous solves */
53: nep->reason = NEP_CONVERGED_ITERATING;
55: /* set default solver type (NEPSetFromOptions was not called) */
56: if (!((PetscObject)nep)->type_name) {
57: NEPSetType(nep,NEPRII);
58: }
59: if (nep->useds && !nep->ds) { NEPGetDS(nep,&nep->ds); }
60: if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
61: if (!((PetscObject)nep->rg)->type_name) {
62: RGSetType(nep->rg,RGINTERVAL);
63: }
64: if (nep->twosided && !nep->hasts) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"This solver does not support computing left eigenvectors (no two-sided variant)");
66: /* set problem dimensions */
67: switch (nep->fui) {
68: case NEP_USER_INTERFACE_CALLBACK:
69: NEPGetFunction(nep,&T,NULL,NULL,NULL);
70: MatGetSize(T,&nep->n,NULL);
71: MatGetLocalSize(T,&nep->nloc,NULL);
72: break;
73: case NEP_USER_INTERFACE_SPLIT:
74: MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->function);
75: MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->jacobian);
76: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->function);
77: PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->jacobian);
78: MatGetSize(nep->A[0],&nep->n,NULL);
79: MatGetLocalSize(nep->A[0],&nep->nloc,NULL);
80: break;
81: }
83: /* set default problem type */
84: if (!nep->problem_type) {
85: NEPSetProblemType(nep,NEP_GENERAL);
86: }
88: /* check consistency of refinement options */
89: if (nep->refine) {
90: if (nep->fui!=NEP_USER_INTERFACE_SPLIT) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Iterative refinement only implemented in split form");
91: if (!nep->scheme) { /* set default scheme */
92: NEPRefineGetKSP(nep,&ksp);
93: KSPGetPC(ksp,&pc);
94: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
95: if (flg) {
96: PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
97: }
98: nep->scheme = flg? NEP_REFINE_SCHEME_MBE: NEP_REFINE_SCHEME_SCHUR;
99: }
100: if (nep->scheme==NEP_REFINE_SCHEME_MBE) {
101: NEPRefineGetKSP(nep,&ksp);
102: KSPGetPC(ksp,&pc);
103: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
104: if (flg) {
105: PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
106: }
107: if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The MBE scheme for refinement requires a direct solver in KSP");
108: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
109: if (size>1) { /* currently selected PC is a factorization */
110: PCFactorGetMatSolverType(pc,&stype);
111: PetscStrcmp(stype,MATSOLVERPETSC,&flg);
112: if (flg) SETERRQ(PetscObjectComm((PetscObject)nep),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");
113: }
114: }
115: if (nep->scheme==NEP_REFINE_SCHEME_SCHUR) {
116: if (nep->npart>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The Schur scheme for refinement does not support subcommunicators");
117: }
118: }
119: /* call specific solver setup */
120: (*nep->ops->setup)(nep);
122: /* by default, compute eigenvalues close to target */
123: /* nep->target should contain the initial guess for the eigenvalue */
124: if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;
126: /* set tolerance if not yet set */
127: if (nep->tol==PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
128: if (nep->refine) {
129: if (nep->rtol==PETSC_DEFAULT) nep->rtol = PetscMax(nep->tol/1000,PETSC_MACHINE_EPSILON);
130: if (nep->rits==PETSC_DEFAULT) nep->rits = (nep->refine==NEP_REFINE_SIMPLE)? 10: 1;
131: }
133: /* fill sorting criterion context */
134: switch (nep->which) {
135: case NEP_LARGEST_MAGNITUDE:
136: nep->sc->comparison = SlepcCompareLargestMagnitude;
137: nep->sc->comparisonctx = NULL;
138: break;
139: case NEP_SMALLEST_MAGNITUDE:
140: nep->sc->comparison = SlepcCompareSmallestMagnitude;
141: nep->sc->comparisonctx = NULL;
142: break;
143: case NEP_LARGEST_REAL:
144: nep->sc->comparison = SlepcCompareLargestReal;
145: nep->sc->comparisonctx = NULL;
146: break;
147: case NEP_SMALLEST_REAL:
148: nep->sc->comparison = SlepcCompareSmallestReal;
149: nep->sc->comparisonctx = NULL;
150: break;
151: case NEP_LARGEST_IMAGINARY:
152: nep->sc->comparison = SlepcCompareLargestImaginary;
153: nep->sc->comparisonctx = NULL;
154: break;
155: case NEP_SMALLEST_IMAGINARY:
156: nep->sc->comparison = SlepcCompareSmallestImaginary;
157: nep->sc->comparisonctx = NULL;
158: break;
159: case NEP_TARGET_MAGNITUDE:
160: nep->sc->comparison = SlepcCompareTargetMagnitude;
161: nep->sc->comparisonctx = &nep->target;
162: break;
163: case NEP_TARGET_REAL:
164: nep->sc->comparison = SlepcCompareTargetReal;
165: nep->sc->comparisonctx = &nep->target;
166: break;
167: case NEP_TARGET_IMAGINARY:
168: #if defined(PETSC_USE_COMPLEX)
169: nep->sc->comparison = SlepcCompareTargetImaginary;
170: nep->sc->comparisonctx = &nep->target;
171: #endif
172: break;
173: case NEP_ALL:
174: nep->sc->comparison = SlepcCompareSmallestReal;
175: nep->sc->comparisonctx = NULL;
176: break;
177: case NEP_WHICH_USER:
178: break;
179: }
181: nep->sc->map = NULL;
182: nep->sc->mapobj = NULL;
184: /* fill sorting criterion for DS */
185: if (nep->useds) {
186: DSGetSlepcSC(nep->ds,&sc);
187: sc->comparison = nep->sc->comparison;
188: sc->comparisonctx = nep->sc->comparisonctx;
189: PetscObjectTypeCompare((PetscObject)nep,NEPNLEIGS,&flg);
190: if (!flg) {
191: sc->map = NULL;
192: sc->mapobj = NULL;
193: }
194: }
195: if (nep->nev > nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");
197: /* process initial vectors */
198: if (nep->nini<0) {
199: k = -nep->nini;
200: if (k>nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The number of initial vectors is larger than ncv");
201: BVInsertVecs(nep->V,0,&k,nep->IS,PETSC_TRUE);
202: SlepcBasisDestroy_Private(&nep->nini,&nep->IS);
203: nep->nini = k;
204: }
205: PetscLogEventEnd(NEP_SetUp,nep,0,0,0);
206: nep->state = NEP_STATE_SETUP;
207: return(0);
208: }
210: /*@C
211: NEPSetInitialSpace - Specify a basis of vectors that constitute the initial
212: space, that is, the subspace from which the solver starts to iterate.
214: Collective on nep
216: Input Parameter:
217: + nep - the nonlinear eigensolver context
218: . n - number of vectors
219: - is - set of basis vectors of the initial space
221: Notes:
222: Some solvers start to iterate on a single vector (initial vector). In that case,
223: the other vectors are ignored.
225: These vectors do not persist from one NEPSolve() call to the other, so the
226: initial space should be set every time.
228: The vectors do not need to be mutually orthonormal, since they are explicitly
229: orthonormalized internally.
231: Common usage of this function is when the user can provide a rough approximation
232: of the wanted eigenspace. Then, convergence may be faster.
234: Level: intermediate
235: @*/
236: PetscErrorCode NEPSetInitialSpace(NEP nep,PetscInt n,Vec is[])237: {
243: if (n<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
244: if (n>0) {
247: }
248: SlepcBasisReference_Private(n,is,&nep->nini,&nep->IS);
249: if (n>0) nep->state = NEP_STATE_INITIAL;
250: return(0);
251: }
253: /*
254: NEPSetDimensions_Default - Set reasonable values for ncv, mpd if not set
255: by the user. This is called at setup.
256: */
257: PetscErrorCode NEPSetDimensions_Default(NEP nep,PetscInt nev,PetscInt *ncv,PetscInt *mpd)258: {
260: if (*ncv) { /* ncv set */
261: if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must be at least nev");
262: } else if (*mpd) { /* mpd set */
263: *ncv = PetscMin(nep->n,nev+(*mpd));
264: } else { /* neither set: defaults depend on nev being small or large */
265: if (nev<500) *ncv = PetscMin(nep->n,PetscMax(2*nev,nev+15));
266: else {
267: *mpd = 500;
268: *ncv = PetscMin(nep->n,nev+(*mpd));
269: }
270: }
271: if (!*mpd) *mpd = *ncv;
272: return(0);
273: }
275: /*@
276: NEPAllocateSolution - Allocate memory storage for common variables such
277: as eigenvalues and eigenvectors.
279: Collective on nep
281: Input Parameters:
282: + nep - eigensolver context
283: - extra - number of additional positions, used for methods that require a
284: working basis slightly larger than ncv
286: Developers Note:
287: This is SLEPC_EXTERN because it may be required by user plugin NEP288: implementations.
290: Level: developer
291: @*/
292: PetscErrorCode NEPAllocateSolution(NEP nep,PetscInt extra)293: {
295: PetscInt oldsize,newc,requested;
296: PetscLogDouble cnt;
297: Mat T;
298: Vec t;
301: requested = nep->ncv + extra;
303: /* oldsize is zero if this is the first time setup is called */
304: BVGetSizes(nep->V,NULL,NULL,&oldsize);
305: newc = PetscMax(0,requested-oldsize);
307: /* allocate space for eigenvalues and friends */
308: if (requested != oldsize || !nep->eigr) {
309: PetscFree4(nep->eigr,nep->eigi,nep->errest,nep->perm);
310: PetscMalloc4(requested,&nep->eigr,requested,&nep->eigi,requested,&nep->errest,requested,&nep->perm);
311: cnt = newc*sizeof(PetscScalar) + newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
312: PetscLogObjectMemory((PetscObject)nep,cnt);
313: }
315: /* allocate V */
316: if (!nep->V) { NEPGetBV(nep,&nep->V); }
317: if (!oldsize) {
318: if (!((PetscObject)(nep->V))->type_name) {
319: BVSetType(nep->V,BVSVEC);
320: }
321: if (nep->fui==NEP_USER_INTERFACE_SPLIT) T = nep->A[0];
322: else {
323: NEPGetFunction(nep,&T,NULL,NULL,NULL);
324: }
325: MatCreateVecsEmpty(T,&t,NULL);
326: BVSetSizesFromVec(nep->V,t,requested);
327: VecDestroy(&t);
328: } else {
329: BVResize(nep->V,requested,PETSC_FALSE);
330: }
332: /* allocate W */
333: if (nep->twosided) {
334: BVDestroy(&nep->W);
335: BVDuplicate(nep->V,&nep->W);
336: }
337: return(0);
338: }