Actual source code: nepsetup.c

slepc-3.16.2 2022-02-01
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2021, 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: {
 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:   }

 65:   /* set problem dimensions */
 66:   switch (nep->fui) {
 67:   case NEP_USER_INTERFACE_CALLBACK:
 68:     NEPGetFunction(nep,&T,NULL,NULL,NULL);
 69:     MatGetSize(T,&nep->n,NULL);
 70:     MatGetLocalSize(T,&nep->nloc,NULL);
 71:     break;
 72:   case NEP_USER_INTERFACE_SPLIT:
 73:     MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->function);
 74:     MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->jacobian);
 75:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->function);
 76:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->jacobian);
 77:     MatGetSize(nep->A[0],&nep->n,NULL);
 78:     MatGetLocalSize(nep->A[0],&nep->nloc,NULL);
 79:     break;
 80:   }

 82:   /* set default problem type */
 83:   if (!nep->problem_type) {
 84:     NEPSetProblemType(nep,NEP_GENERAL);
 85:   }

 87:   /* check consistency of refinement options */
 88:   if (nep->refine) {
 89:     if (nep->fui!=NEP_USER_INTERFACE_SPLIT) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Iterative refinement only implemented in split form");
 90:     if (!nep->scheme) {  /* set default scheme */
 91:       NEPRefineGetKSP(nep,&ksp);
 92:       KSPGetPC(ksp,&pc);
 93:       PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
 94:       if (flg) {
 95:         PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
 96:       }
 97:       nep->scheme = flg? NEP_REFINE_SCHEME_MBE: NEP_REFINE_SCHEME_SCHUR;
 98:     }
 99:     if (nep->scheme==NEP_REFINE_SCHEME_MBE) {
100:       NEPRefineGetKSP(nep,&ksp);
101:       KSPGetPC(ksp,&pc);
102:       PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
103:       if (flg) {
104:         PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
105:       }
106:       if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The MBE scheme for refinement requires a direct solver in KSP");
107:       MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
108:       if (size>1) {   /* currently selected PC is a factorization */
109:         PCFactorGetMatSolverType(pc,&stype);
110:         PetscStrcmp(stype,MATSOLVERPETSC,&flg);
111:         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");
112:       }
113:     }
114:     if (nep->scheme==NEP_REFINE_SCHEME_SCHUR) {
115:       if (nep->npart>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The Schur scheme for refinement does not support subcommunicators");
116:     }
117:   }
118:   /* call specific solver setup */
119:   (*nep->ops->setup)(nep);

121:   /* set tolerance if not yet set */
122:   if (nep->tol==PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
123:   if (nep->refine) {
124:     if (nep->rtol==PETSC_DEFAULT) nep->rtol = PetscMax(nep->tol/1000,PETSC_MACHINE_EPSILON);
125:     if (nep->rits==PETSC_DEFAULT) nep->rits = (nep->refine==NEP_REFINE_SIMPLE)? 10: 1;
126:   }

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

176:   nep->sc->map    = NULL;
177:   nep->sc->mapobj = NULL;

179:   /* fill sorting criterion for DS */
180:   if (nep->useds) {
181:     DSGetSlepcSC(nep->ds,&sc);
182:     sc->comparison    = nep->sc->comparison;
183:     sc->comparisonctx = nep->sc->comparisonctx;
184:     PetscObjectTypeCompare((PetscObject)nep,NEPNLEIGS,&flg);
185:     if (!flg) {
186:       sc->map    = NULL;
187:       sc->mapobj = NULL;
188:     }
189:   }
190:   if (nep->nev > nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");

192:   /* process initial vectors */
193:   if (nep->nini<0) {
194:     k = -nep->nini;
195:     if (k>nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER_INPUT,"The number of initial vectors is larger than ncv");
196:     BVInsertVecs(nep->V,0,&k,nep->IS,PETSC_TRUE);
197:     SlepcBasisDestroy_Private(&nep->nini,&nep->IS);
198:     nep->nini = k;
199:   }
200:   PetscLogEventEnd(NEP_SetUp,nep,0,0,0);
201:   nep->state = NEP_STATE_SETUP;
202:   return(0);
203: }

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

209:    Collective on nep

211:    Input Parameters:
212: +  nep   - the nonlinear eigensolver context
213: .  n     - number of vectors
214: -  is    - set of basis vectors of the initial space

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

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

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

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

229:    Level: intermediate
230: @*/
231: PetscErrorCode NEPSetInitialSpace(NEP nep,PetscInt n,Vec is[])
232: {

238:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
239:   if (n>0) {
242:   }
243:   SlepcBasisReference_Private(n,is,&nep->nini,&nep->IS);
244:   if (n>0) nep->state = NEP_STATE_INITIAL;
245:   return(0);
246: }

248: /*
249:   NEPSetDimensions_Default - Set reasonable values for ncv, mpd if not set
250:   by the user. This is called at setup.
251:  */
252: PetscErrorCode NEPSetDimensions_Default(NEP nep,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
253: {
255:   if (*ncv!=PETSC_DEFAULT) { /* ncv set */
256:     if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_USER_INPUT,"The value of ncv must be at least nev");
257:   } else if (*mpd!=PETSC_DEFAULT) { /* mpd set */
258:     *ncv = PetscMin(nep->n,nev+(*mpd));
259:   } else { /* neither set: defaults depend on nev being small or large */
260:     if (nev<500) *ncv = PetscMin(nep->n,PetscMax(2*nev,nev+15));
261:     else {
262:       *mpd = 500;
263:       *ncv = PetscMin(nep->n,nev+(*mpd));
264:     }
265:   }
266:   if (*mpd==PETSC_DEFAULT) *mpd = *ncv;
267:   return(0);
268: }

270: /*@
271:    NEPAllocateSolution - Allocate memory storage for common variables such
272:    as eigenvalues and eigenvectors.

274:    Collective on nep

276:    Input Parameters:
277: +  nep   - eigensolver context
278: -  extra - number of additional positions, used for methods that require a
279:            working basis slightly larger than ncv

281:    Developers Note:
282:    This is SLEPC_EXTERN because it may be required by user plugin NEP
283:    implementations.

285:    Level: developer
286: @*/
287: PetscErrorCode NEPAllocateSolution(NEP nep,PetscInt extra)
288: {
290:   PetscInt       oldsize,newc,requested;
291:   PetscLogDouble cnt;
292:   PetscRandom    rand;
293:   Mat            T;
294:   Vec            t;

297:   requested = nep->ncv + extra;

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

303:   /* allocate space for eigenvalues and friends */
304:   if (requested != oldsize || !nep->eigr) {
305:     PetscFree4(nep->eigr,nep->eigi,nep->errest,nep->perm);
306:     PetscMalloc4(requested,&nep->eigr,requested,&nep->eigi,requested,&nep->errest,requested,&nep->perm);
307:     cnt = newc*sizeof(PetscScalar) + newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
308:     PetscLogObjectMemory((PetscObject)nep,cnt);
309:   }

311:   /* allocate V */
312:   if (!nep->V) { NEPGetBV(nep,&nep->V); }
313:   if (!oldsize) {
314:     if (!((PetscObject)(nep->V))->type_name) {
315:       BVSetType(nep->V,BVSVEC);
316:     }
317:     if (nep->fui==NEP_USER_INTERFACE_SPLIT) T = nep->A[0];
318:     else {
319:       NEPGetFunction(nep,&T,NULL,NULL,NULL);
320:     }
321:     MatCreateVecsEmpty(T,&t,NULL);
322:     BVSetSizesFromVec(nep->V,t,requested);
323:     VecDestroy(&t);
324:   } else {
325:     BVResize(nep->V,requested,PETSC_FALSE);
326:   }

328:   /* allocate W */
329:   if (nep->twosided) {
330:     BVGetRandomContext(nep->V,&rand);  /* make sure the random context is available when duplicating */
331:     BVDestroy(&nep->W);
332:     BVDuplicate(nep->V,&nep->W);
333:   }
334:   return(0);
335: }