Actual source code: nepsetup.c
slepc-3.18.0 2022-10-01
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: NEP routines related to problem setup
12: */
14: #include <slepc/private/nepimpl.h>
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: {
36: PetscInt k;
37: SlepcSC sc;
38: Mat T;
39: PetscBool flg;
40: KSP ksp;
41: PC pc;
42: PetscMPIInt size;
43: MatSolverType stype;
46: NEPCheckProblem(nep,1);
47: if (nep->state) return 0;
48: PetscLogEventBegin(NEP_SetUp,nep,0,0,0);
50: /* reset the convergence flag from the previous solves */
51: nep->reason = NEP_CONVERGED_ITERATING;
53: /* set default solver type (NEPSetFromOptions was not called) */
54: if (!((PetscObject)nep)->type_name) NEPSetType(nep,NEPRII);
55: if (nep->useds && !nep->ds) NEPGetDS(nep,&nep->ds);
56: if (!nep->rg) NEPGetRG(nep,&nep->rg);
57: if (!((PetscObject)nep->rg)->type_name) RGSetType(nep->rg,RGINTERVAL);
59: /* set problem dimensions */
60: switch (nep->fui) {
61: case NEP_USER_INTERFACE_CALLBACK:
62: NEPGetFunction(nep,&T,NULL,NULL,NULL);
63: MatGetSize(T,&nep->n,NULL);
64: MatGetLocalSize(T,&nep->nloc,NULL);
65: break;
66: case NEP_USER_INTERFACE_SPLIT:
67: MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->function);
68: if (nep->P) MatDuplicate(nep->P[0],MAT_DO_NOT_COPY_VALUES,&nep->function_pre);
69: MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->jacobian);
70: MatGetSize(nep->A[0],&nep->n,NULL);
71: MatGetLocalSize(nep->A[0],&nep->nloc,NULL);
72: break;
73: }
75: /* set default problem type */
76: if (!nep->problem_type) NEPSetProblemType(nep,NEP_GENERAL);
78: /* check consistency of refinement options */
79: if (nep->refine) {
81: if (!nep->scheme) { /* set default scheme */
82: NEPRefineGetKSP(nep,&ksp);
83: KSPGetPC(ksp,&pc);
84: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
85: if (flg) PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
86: nep->scheme = flg? NEP_REFINE_SCHEME_MBE: NEP_REFINE_SCHEME_SCHUR;
87: }
88: if (nep->scheme==NEP_REFINE_SCHEME_MBE) {
89: NEPRefineGetKSP(nep,&ksp);
90: KSPGetPC(ksp,&pc);
91: PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
92: if (flg) PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
94: MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
95: if (size>1) { /* currently selected PC is a factorization */
96: PCFactorGetMatSolverType(pc,&stype);
97: PetscStrcmp(stype,MATSOLVERPETSC,&flg);
99: }
100: }
101: if (nep->scheme==NEP_REFINE_SCHEME_SCHUR) {
103: }
104: }
105: /* call specific solver setup */
106: PetscUseTypeMethod(nep,setup);
108: /* set tolerance if not yet set */
109: if (nep->tol==PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
110: if (nep->refine) {
111: if (nep->rtol==PETSC_DEFAULT) nep->rtol = PetscMax(nep->tol/1000,PETSC_MACHINE_EPSILON);
112: if (nep->rits==PETSC_DEFAULT) nep->rits = (nep->refine==NEP_REFINE_SIMPLE)? 10: 1;
113: }
115: /* fill sorting criterion context */
116: switch (nep->which) {
117: case NEP_LARGEST_MAGNITUDE:
118: nep->sc->comparison = SlepcCompareLargestMagnitude;
119: nep->sc->comparisonctx = NULL;
120: break;
121: case NEP_SMALLEST_MAGNITUDE:
122: nep->sc->comparison = SlepcCompareSmallestMagnitude;
123: nep->sc->comparisonctx = NULL;
124: break;
125: case NEP_LARGEST_REAL:
126: nep->sc->comparison = SlepcCompareLargestReal;
127: nep->sc->comparisonctx = NULL;
128: break;
129: case NEP_SMALLEST_REAL:
130: nep->sc->comparison = SlepcCompareSmallestReal;
131: nep->sc->comparisonctx = NULL;
132: break;
133: case NEP_LARGEST_IMAGINARY:
134: nep->sc->comparison = SlepcCompareLargestImaginary;
135: nep->sc->comparisonctx = NULL;
136: break;
137: case NEP_SMALLEST_IMAGINARY:
138: nep->sc->comparison = SlepcCompareSmallestImaginary;
139: nep->sc->comparisonctx = NULL;
140: break;
141: case NEP_TARGET_MAGNITUDE:
142: nep->sc->comparison = SlepcCompareTargetMagnitude;
143: nep->sc->comparisonctx = &nep->target;
144: break;
145: case NEP_TARGET_REAL:
146: nep->sc->comparison = SlepcCompareTargetReal;
147: nep->sc->comparisonctx = &nep->target;
148: break;
149: case NEP_TARGET_IMAGINARY:
150: #if defined(PETSC_USE_COMPLEX)
151: nep->sc->comparison = SlepcCompareTargetImaginary;
152: nep->sc->comparisonctx = &nep->target;
153: #endif
154: break;
155: case NEP_ALL:
156: nep->sc->comparison = SlepcCompareSmallestReal;
157: nep->sc->comparisonctx = NULL;
158: break;
159: case NEP_WHICH_USER:
160: break;
161: }
163: nep->sc->map = NULL;
164: nep->sc->mapobj = NULL;
166: /* fill sorting criterion for DS */
167: if (nep->useds) {
168: DSGetSlepcSC(nep->ds,&sc);
169: sc->comparison = nep->sc->comparison;
170: sc->comparisonctx = nep->sc->comparisonctx;
171: PetscObjectTypeCompare((PetscObject)nep,NEPNLEIGS,&flg);
172: if (!flg) {
173: sc->map = NULL;
174: sc->mapobj = NULL;
175: }
176: }
179: /* process initial vectors */
180: if (nep->nini<0) {
181: k = -nep->nini;
183: BVInsertVecs(nep->V,0,&k,nep->IS,PETSC_TRUE);
184: SlepcBasisDestroy_Private(&nep->nini,&nep->IS);
185: nep->nini = k;
186: }
187: PetscLogEventEnd(NEP_SetUp,nep,0,0,0);
188: nep->state = NEP_STATE_SETUP;
189: return 0;
190: }
192: /*@C
193: NEPSetInitialSpace - Specify a basis of vectors that constitute the initial
194: space, that is, the subspace from which the solver starts to iterate.
196: Collective on nep
198: Input Parameters:
199: + nep - the nonlinear eigensolver context
200: . n - number of vectors
201: - is - set of basis vectors of the initial space
203: Notes:
204: Some solvers start to iterate on a single vector (initial vector). In that case,
205: the other vectors are ignored.
207: These vectors do not persist from one NEPSolve() call to the other, so the
208: initial space should be set every time.
210: The vectors do not need to be mutually orthonormal, since they are explicitly
211: orthonormalized internally.
213: Common usage of this function is when the user can provide a rough approximation
214: of the wanted eigenspace. Then, convergence may be faster.
216: Level: intermediate
218: .seealso: NEPSetUp()
219: @*/
220: PetscErrorCode NEPSetInitialSpace(NEP nep,PetscInt n,Vec is[])
221: {
225: if (n>0) {
228: }
229: SlepcBasisReference_Private(n,is,&nep->nini,&nep->IS);
230: if (n>0) nep->state = NEP_STATE_INITIAL;
231: return 0;
232: }
234: /*
235: NEPSetDimensions_Default - Set reasonable values for ncv, mpd if not set
236: by the user. This is called at setup.
237: */
238: PetscErrorCode NEPSetDimensions_Default(NEP nep,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
239: {
240: if (*ncv!=PETSC_DEFAULT) { /* ncv set */
242: } else if (*mpd!=PETSC_DEFAULT) { /* mpd set */
243: *ncv = PetscMin(nep->n,nev+(*mpd));
244: } else { /* neither set: defaults depend on nev being small or large */
245: if (nev<500) *ncv = PetscMin(nep->n,PetscMax(2*nev,nev+15));
246: else {
247: *mpd = 500;
248: *ncv = PetscMin(nep->n,nev+(*mpd));
249: }
250: }
251: if (*mpd==PETSC_DEFAULT) *mpd = *ncv;
252: return 0;
253: }
255: /*@
256: NEPAllocateSolution - Allocate memory storage for common variables such
257: as eigenvalues and eigenvectors.
259: Collective on nep
261: Input Parameters:
262: + nep - eigensolver context
263: - extra - number of additional positions, used for methods that require a
264: working basis slightly larger than ncv
266: Developer Notes:
267: This is SLEPC_EXTERN because it may be required by user plugin NEP
268: implementations.
270: Level: developer
272: .seealso: PEPSetUp()
273: @*/
274: PetscErrorCode NEPAllocateSolution(NEP nep,PetscInt extra)
275: {
276: PetscInt oldsize,requested;
277: PetscRandom rand;
278: Mat T;
279: Vec t;
281: requested = nep->ncv + extra;
283: /* oldsize is zero if this is the first time setup is called */
284: BVGetSizes(nep->V,NULL,NULL,&oldsize);
286: /* allocate space for eigenvalues and friends */
287: if (requested != oldsize || !nep->eigr) {
288: PetscFree4(nep->eigr,nep->eigi,nep->errest,nep->perm);
289: PetscMalloc4(requested,&nep->eigr,requested,&nep->eigi,requested,&nep->errest,requested,&nep->perm);
290: }
292: /* allocate V */
293: if (!nep->V) NEPGetBV(nep,&nep->V);
294: if (!oldsize) {
295: if (!((PetscObject)(nep->V))->type_name) BVSetType(nep->V,BVSVEC);
296: if (nep->fui==NEP_USER_INTERFACE_SPLIT) T = nep->A[0];
297: else NEPGetFunction(nep,&T,NULL,NULL,NULL);
298: MatCreateVecsEmpty(T,&t,NULL);
299: BVSetSizesFromVec(nep->V,t,requested);
300: VecDestroy(&t);
301: } else BVResize(nep->V,requested,PETSC_FALSE);
303: /* allocate W */
304: if (nep->twosided) {
305: BVGetRandomContext(nep->V,&rand); /* make sure the random context is available when duplicating */
306: BVDestroy(&nep->W);
307: BVDuplicate(nep->V,&nep->W);
308: }
309: return 0;
310: }