Actual source code: nepsetup.c

slepc-3.10.1 2018-10-23
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2018, 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->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:   case NEP_USER_INTERFACE_DERIVATIVES:
 81:     NEPGetDerivatives(nep,&T,NULL,NULL);
 82:     MatGetSize(T,&nep->n,NULL);
 83:     MatGetLocalSize(T,&nep->nloc,NULL);
 84:     break;
 85:   }

 87:   /* set default problem type */
 88:   if (!nep->problem_type) {
 89:     NEPSetProblemType(nep,NEP_GENERAL);
 90:   }

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

126:   /* by default, compute eigenvalues close to target */
127:   /* nep->target should contain the initial guess for the eigenvalue */
128:   if (!nep->which) nep->which = NEP_TARGET_MAGNITUDE;

130:   /* set tolerance if not yet set */
131:   if (nep->tol==PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
132:   if (nep->refine) {
133:     if (nep->rtol==PETSC_DEFAULT) nep->rtol = PetscMax(nep->tol/1000,PETSC_MACHINE_EPSILON);
134:     if (nep->rits==PETSC_DEFAULT) nep->rits = (nep->refine==NEP_REFINE_SIMPLE)? 10: 1;
135:   }

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

185:   nep->sc->map    = NULL;
186:   nep->sc->mapobj = NULL;

188:   /* fill sorting criterion for DS */
189:   DSGetSlepcSC(nep->ds,&sc);
190:   sc->comparison    = nep->sc->comparison;
191:   sc->comparisonctx = nep->sc->comparisonctx;
192:   PetscObjectTypeCompare((PetscObject)nep,NEPNLEIGS,&flg);
193:   if (!flg) {
194:     sc->map    = NULL;
195:     sc->mapobj = NULL;
196:   }
197:   if (nep->nev > nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");

199:   /* process initial vectors */
200:   if (nep->nini<0) {
201:     k = -nep->nini;
202:     if (k>nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The number of initial vectors is larger than ncv");
203:     BVInsertVecs(nep->V,0,&k,nep->IS,PETSC_TRUE);
204:     SlepcBasisDestroy_Private(&nep->nini,&nep->IS);
205:     nep->nini = k;
206:   }
207:   PetscLogEventEnd(NEP_SetUp,nep,0,0,0);
208:   nep->state = NEP_STATE_SETUP;
209:   return(0);
210: }

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

216:    Collective on NEP and Vec

218:    Input Parameter:
219: +  nep   - the nonlinear eigensolver context
220: .  n     - number of vectors
221: -  is    - set of basis vectors of the initial space

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

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

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

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

236:    Level: intermediate
237: @*/
238: PetscErrorCode NEPSetInitialSpace(NEP nep,PetscInt n,Vec *is)
239: {

245:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
246:   SlepcBasisReference_Private(n,is,&nep->nini,&nep->IS);
247:   if (n>0) nep->state = NEP_STATE_INITIAL;
248:   return(0);
249: }

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

273: /*@
274:    NEPAllocateSolution - Allocate memory storage for common variables such
275:    as eigenvalues and eigenvectors.

277:    Collective on NEP

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

284:    Developers Note:
285:    This is PETSC_EXTERN because it may be required by user plugin NEP
286:    implementations.

288:    Level: developer
289: @*/
290: PetscErrorCode NEPAllocateSolution(NEP nep,PetscInt extra)
291: {
293:   PetscInt       oldsize,newc,requested;
294:   PetscLogDouble cnt;
295:   Mat            T;
296:   Vec            t;

299:   requested = nep->ncv + extra;

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

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

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