Actual source code: nepsetup.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  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) PetscFunctionReturn(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:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->function);
 69:     if (nep->P) {
 70:       MatDuplicate(nep->P[0],MAT_DO_NOT_COPY_VALUES,&nep->function_pre);
 71:       PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->function_pre);
 72:     }
 73:     MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->jacobian);
 74:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->jacobian);
 75:     MatGetSize(nep->A[0],&nep->n,NULL);
 76:     MatGetLocalSize(nep->A[0],&nep->nloc,NULL);
 77:     break;
 78:   }

 80:   /* set default problem type */
 81:   if (!nep->problem_type) NEPSetProblemType(nep,NEP_GENERAL);

 83:   /* check consistency of refinement options */
 84:   if (nep->refine) {
 86:     if (!nep->scheme) {  /* set default scheme */
 87:       NEPRefineGetKSP(nep,&ksp);
 88:       KSPGetPC(ksp,&pc);
 89:       PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
 90:       if (flg) PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
 91:       nep->scheme = flg? NEP_REFINE_SCHEME_MBE: NEP_REFINE_SCHEME_SCHUR;
 92:     }
 93:     if (nep->scheme==NEP_REFINE_SCHEME_MBE) {
 94:       NEPRefineGetKSP(nep,&ksp);
 95:       KSPGetPC(ksp,&pc);
 96:       PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
 97:       if (flg) PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
 99:       MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
100:       if (size>1) {   /* currently selected PC is a factorization */
101:         PCFactorGetMatSolverType(pc,&stype);
102:         PetscStrcmp(stype,MATSOLVERPETSC,&flg);
104:       }
105:     }
106:     if (nep->scheme==NEP_REFINE_SCHEME_SCHUR) {
108:     }
109:   }
110:   /* call specific solver setup */
111:   (*nep->ops->setup)(nep);

113:   /* set tolerance if not yet set */
114:   if (nep->tol==PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
115:   if (nep->refine) {
116:     if (nep->rtol==PETSC_DEFAULT) nep->rtol = PetscMax(nep->tol/1000,PETSC_MACHINE_EPSILON);
117:     if (nep->rits==PETSC_DEFAULT) nep->rits = (nep->refine==NEP_REFINE_SIMPLE)? 10: 1;
118:   }

120:   /* fill sorting criterion context */
121:   switch (nep->which) {
122:     case NEP_LARGEST_MAGNITUDE:
123:       nep->sc->comparison    = SlepcCompareLargestMagnitude;
124:       nep->sc->comparisonctx = NULL;
125:       break;
126:     case NEP_SMALLEST_MAGNITUDE:
127:       nep->sc->comparison    = SlepcCompareSmallestMagnitude;
128:       nep->sc->comparisonctx = NULL;
129:       break;
130:     case NEP_LARGEST_REAL:
131:       nep->sc->comparison    = SlepcCompareLargestReal;
132:       nep->sc->comparisonctx = NULL;
133:       break;
134:     case NEP_SMALLEST_REAL:
135:       nep->sc->comparison    = SlepcCompareSmallestReal;
136:       nep->sc->comparisonctx = NULL;
137:       break;
138:     case NEP_LARGEST_IMAGINARY:
139:       nep->sc->comparison    = SlepcCompareLargestImaginary;
140:       nep->sc->comparisonctx = NULL;
141:       break;
142:     case NEP_SMALLEST_IMAGINARY:
143:       nep->sc->comparison    = SlepcCompareSmallestImaginary;
144:       nep->sc->comparisonctx = NULL;
145:       break;
146:     case NEP_TARGET_MAGNITUDE:
147:       nep->sc->comparison    = SlepcCompareTargetMagnitude;
148:       nep->sc->comparisonctx = &nep->target;
149:       break;
150:     case NEP_TARGET_REAL:
151:       nep->sc->comparison    = SlepcCompareTargetReal;
152:       nep->sc->comparisonctx = &nep->target;
153:       break;
154:     case NEP_TARGET_IMAGINARY:
155: #if defined(PETSC_USE_COMPLEX)
156:       nep->sc->comparison    = SlepcCompareTargetImaginary;
157:       nep->sc->comparisonctx = &nep->target;
158: #endif
159:       break;
160:     case NEP_ALL:
161:       nep->sc->comparison    = SlepcCompareSmallestReal;
162:       nep->sc->comparisonctx = NULL;
163:       break;
164:     case NEP_WHICH_USER:
165:       break;
166:   }

168:   nep->sc->map    = NULL;
169:   nep->sc->mapobj = NULL;

171:   /* fill sorting criterion for DS */
172:   if (nep->useds) {
173:     DSGetSlepcSC(nep->ds,&sc);
174:     sc->comparison    = nep->sc->comparison;
175:     sc->comparisonctx = nep->sc->comparisonctx;
176:     PetscObjectTypeCompare((PetscObject)nep,NEPNLEIGS,&flg);
177:     if (!flg) {
178:       sc->map    = NULL;
179:       sc->mapobj = NULL;
180:     }
181:   }

184:   /* process initial vectors */
185:   if (nep->nini<0) {
186:     k = -nep->nini;
188:     BVInsertVecs(nep->V,0,&k,nep->IS,PETSC_TRUE);
189:     SlepcBasisDestroy_Private(&nep->nini,&nep->IS);
190:     nep->nini = k;
191:   }
192:   PetscLogEventEnd(NEP_SetUp,nep,0,0,0);
193:   nep->state = NEP_STATE_SETUP;
194:   PetscFunctionReturn(0);
195: }

197: /*@C
198:    NEPSetInitialSpace - Specify a basis of vectors that constitute the initial
199:    space, that is, the subspace from which the solver starts to iterate.

201:    Collective on nep

203:    Input Parameters:
204: +  nep   - the nonlinear eigensolver context
205: .  n     - number of vectors
206: -  is    - set of basis vectors of the initial space

208:    Notes:
209:    Some solvers start to iterate on a single vector (initial vector). In that case,
210:    the other vectors are ignored.

212:    These vectors do not persist from one NEPSolve() call to the other, so the
213:    initial space should be set every time.

215:    The vectors do not need to be mutually orthonormal, since they are explicitly
216:    orthonormalized internally.

218:    Common usage of this function is when the user can provide a rough approximation
219:    of the wanted eigenspace. Then, convergence may be faster.

221:    Level: intermediate

223: .seealso: NEPSetUp()
224: @*/
225: PetscErrorCode NEPSetInitialSpace(NEP nep,PetscInt n,Vec is[])
226: {
230:   if (n>0) {
233:   }
234:   SlepcBasisReference_Private(n,is,&nep->nini,&nep->IS);
235:   if (n>0) nep->state = NEP_STATE_INITIAL;
236:   PetscFunctionReturn(0);
237: }

239: /*
240:   NEPSetDimensions_Default - Set reasonable values for ncv, mpd if not set
241:   by the user. This is called at setup.
242:  */
243: PetscErrorCode NEPSetDimensions_Default(NEP nep,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
244: {
245:   if (*ncv!=PETSC_DEFAULT) { /* ncv set */
247:   } else if (*mpd!=PETSC_DEFAULT) { /* mpd set */
248:     *ncv = PetscMin(nep->n,nev+(*mpd));
249:   } else { /* neither set: defaults depend on nev being small or large */
250:     if (nev<500) *ncv = PetscMin(nep->n,PetscMax(2*nev,nev+15));
251:     else {
252:       *mpd = 500;
253:       *ncv = PetscMin(nep->n,nev+(*mpd));
254:     }
255:   }
256:   if (*mpd==PETSC_DEFAULT) *mpd = *ncv;
257:   PetscFunctionReturn(0);
258: }

260: /*@
261:    NEPAllocateSolution - Allocate memory storage for common variables such
262:    as eigenvalues and eigenvectors.

264:    Collective on nep

266:    Input Parameters:
267: +  nep   - eigensolver context
268: -  extra - number of additional positions, used for methods that require a
269:            working basis slightly larger than ncv

271:    Developer Notes:
272:    This is SLEPC_EXTERN because it may be required by user plugin NEP
273:    implementations.

275:    Level: developer

277: .seealso: PEPSetUp()
278: @*/
279: PetscErrorCode NEPAllocateSolution(NEP nep,PetscInt extra)
280: {
281:   PetscInt       oldsize,newc,requested;
282:   PetscLogDouble cnt;
283:   PetscRandom    rand;
284:   Mat            T;
285:   Vec            t;

287:   requested = nep->ncv + extra;

289:   /* oldsize is zero if this is the first time setup is called */
290:   BVGetSizes(nep->V,NULL,NULL,&oldsize);
291:   newc = PetscMax(0,requested-oldsize);

293:   /* allocate space for eigenvalues and friends */
294:   if (requested != oldsize || !nep->eigr) {
295:     PetscFree4(nep->eigr,nep->eigi,nep->errest,nep->perm);
296:     PetscMalloc4(requested,&nep->eigr,requested,&nep->eigi,requested,&nep->errest,requested,&nep->perm);
297:     cnt = newc*sizeof(PetscScalar) + newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
298:     PetscLogObjectMemory((PetscObject)nep,cnt);
299:   }

301:   /* allocate V */
302:   if (!nep->V) NEPGetBV(nep,&nep->V);
303:   if (!oldsize) {
304:     if (!((PetscObject)(nep->V))->type_name) BVSetType(nep->V,BVSVEC);
305:     if (nep->fui==NEP_USER_INTERFACE_SPLIT) T = nep->A[0];
306:     else NEPGetFunction(nep,&T,NULL,NULL,NULL);
307:     MatCreateVecsEmpty(T,&t,NULL);
308:     BVSetSizesFromVec(nep->V,t,requested);
309:     VecDestroy(&t);
310:   } else BVResize(nep->V,requested,PETSC_FALSE);

312:   /* allocate W */
313:   if (nep->twosided) {
314:     BVGetRandomContext(nep->V,&rand);  /* make sure the random context is available when duplicating */
315:     BVDestroy(&nep->W);
316:     BVDuplicate(nep->V,&nep->W);
317:   }
318:   PetscFunctionReturn(0);
319: }