Actual source code: nepsetup.c

slepc-3.7.1 2016-05-27
Report Typos and Errors
  1: /*
  2:       NEP routines related to problem setup.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2016, Universitat Politecnica de Valencia, Spain

  8:    This file is part of SLEPc.

 10:    SLEPc is free software: you can redistribute it and/or modify it under  the
 11:    terms of version 3 of the GNU Lesser General Public License as published by
 12:    the Free Software Foundation.

 14:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 15:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 16:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 17:    more details.

 19:    You  should have received a copy of the GNU Lesser General  Public  License
 20:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 21:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 22: */

 24: #include <slepc/private/nepimpl.h>       /*I "slepcnep.h" I*/

 28: /*@
 29:    NEPSetUp - Sets up all the internal data structures necessary for the
 30:    execution of the NEP solver.

 32:    Collective on NEP

 34:    Input Parameter:
 35: .  nep   - solver context

 37:    Notes:
 38:    This function need not be called explicitly in most cases, since NEPSolve()
 39:    calls it. It can be useful when one wants to measure the set-up time
 40:    separately from the solve time.

 42:    Level: developer

 44: .seealso: NEPCreate(), NEPSolve(), NEPDestroy()
 45: @*/
 46: PetscErrorCode NEPSetUp(NEP nep)
 47: {
 49:   PetscInt       k;
 50:   SlepcSC        sc;
 51:   Mat            T;
 52:   PetscBool      flg;
 53:   KSP            ksp;
 54:   PC             pc;
 55:   PetscMPIInt    size;
 56:   const MatSolverPackage stype;

 60:   NEPCheckProblem(nep,1);
 61:   if (nep->state) return(0);
 62:   PetscLogEventBegin(NEP_SetUp,nep,0,0,0);

 64:   /* reset the convergence flag from the previous solves */
 65:   nep->reason = NEP_CONVERGED_ITERATING;

 67:   /* set default solver type (NEPSetFromOptions was not called) */
 68:   if (!((PetscObject)nep)->type_name) {
 69:     NEPSetType(nep,NEPRII);
 70:   }
 71:   if (!nep->ds) { NEPGetDS(nep,&nep->ds); }
 72:   DSReset(nep->ds);
 73:   if (!nep->rg) { NEPGetRG(nep,&nep->rg); }
 74:   if (!((PetscObject)nep->rg)->type_name) {
 75:     RGSetType(nep->rg,RGINTERVAL);
 76:   }

 78:   /* set problem dimensions */
 79:   switch (nep->fui) {
 80:   case NEP_USER_INTERFACE_CALLBACK:
 81:     NEPGetFunction(nep,&T,NULL,NULL,NULL);
 82:     MatGetSize(T,&nep->n,NULL);
 83:     MatGetLocalSize(T,&nep->nloc,NULL);
 84:     break;
 85:   case NEP_USER_INTERFACE_SPLIT:
 86:     MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->function);
 87:     MatDuplicate(nep->A[0],MAT_DO_NOT_COPY_VALUES,&nep->jacobian);
 88:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->function);
 89:     PetscLogObjectParent((PetscObject)nep,(PetscObject)nep->jacobian);
 90:     MatGetSize(nep->A[0],&nep->n,NULL);
 91:     MatGetLocalSize(nep->A[0],&nep->nloc,NULL);
 92:     break;
 93:   case NEP_USER_INTERFACE_DERIVATIVES:
 94:     NEPGetDerivatives(nep,&T,NULL,NULL);
 95:     MatGetSize(T,&nep->n,NULL);
 96:     MatGetLocalSize(T,&nep->nloc,NULL);
 97:     break;
 98:   }

100:   /* check consistency of refinement options */
101:   if (nep->refine) {
102:     if (nep->fui!=NEP_USER_INTERFACE_SPLIT) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"Iterative refinement only implemented in split form");
103:     if (!nep->scheme) {  /* set default scheme */
104:       NEPRefineGetKSP(nep,&ksp);
105:       KSPGetPC(ksp,&pc);
106:       PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
107:       if (flg) {
108:         PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
109:       }
110:       nep->scheme = flg? NEP_REFINE_SCHEME_MBE: NEP_REFINE_SCHEME_SCHUR;
111:     }
112:     if (nep->scheme==NEP_REFINE_SCHEME_MBE) {
113:       NEPRefineGetKSP(nep,&ksp);
114:       KSPGetPC(ksp,&pc);
115:       PetscObjectTypeCompare((PetscObject)ksp,KSPPREONLY,&flg);
116:       if (flg) {
117:         PetscObjectTypeCompareAny((PetscObject)pc,&flg,PCLU,PCCHOLESKY,"");
118:       }
119:       if (!flg) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The MBE scheme for refinement requires a direct solver in KSP");
120:       MPI_Comm_size(PetscObjectComm((PetscObject)pc),&size);
121:       if (size>1) {   /* currently selected PC is a factorization */
122:         PCFactorGetMatSolverPackage(pc,&stype);
123:         PetscStrcmp(stype,MATSOLVERPETSC,&flg);
124:         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");
125:       }
126:     }
127:     if (nep->scheme==NEP_REFINE_SCHEME_SCHUR) {
128:       if (nep->nt<=2) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The Schur scheme for refinement is only available for problems defined with more than 2 matrices");
129:       if (nep->npart>1) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_SUP,"The Schur scheme for refinement does not support subcommunicators"); 
130:     }
131:   }
132:   /* call specific solver setup */
133:   (*nep->ops->setup)(nep);

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

139:   /* set tolerance if not yet set */
140:   if (nep->tol==PETSC_DEFAULT) nep->tol = SLEPC_DEFAULT_TOL;
141:   if (nep->refine) {
142:     if (nep->rtol==PETSC_DEFAULT) nep->rtol = PetscMax(nep->tol/1000,PETSC_MACHINE_EPSILON);
143:     if (nep->rits==PETSC_DEFAULT) nep->rits = (nep->refine==NEP_REFINE_SIMPLE)? 10: 1;
144:   }

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

192:   nep->sc->map    = NULL;
193:   nep->sc->mapobj = NULL;

195:   /* fill sorting criterion for DS */
196:   DSGetSlepcSC(nep->ds,&sc);
197:   sc->comparison    = nep->sc->comparison;
198:   sc->comparisonctx = nep->sc->comparisonctx;
199:   PetscObjectTypeCompare((PetscObject)nep,NEPNLEIGS,&flg);
200:   if (!flg) {
201:     sc->map    = NULL;
202:     sc->mapobj = NULL;
203:   }
204:   if (nep->ncv > nep->n) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"ncv must be the problem size at most");
205:   if (nep->nev > nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");

207:   /* process initial vectors */
208:   if (nep->nini<0) {
209:     k = -nep->nini;
210:     if (k>nep->ncv) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The number of initial vectors is larger than ncv");
211:     BVInsertVecs(nep->V,0,&k,nep->IS,PETSC_TRUE);
212:     SlepcBasisDestroy_Private(&nep->nini,&nep->IS);
213:     nep->nini = k;
214:   }
215:   PetscLogEventEnd(NEP_SetUp,nep,0,0,0);
216:   nep->state = NEP_STATE_SETUP;
217:   return(0);
218: }

222: /*@
223:    NEPSetInitialSpace - Specify a basis of vectors that constitute the initial
224:    space, that is, the subspace from which the solver starts to iterate.

226:    Collective on NEP and Vec

228:    Input Parameter:
229: +  nep   - the nonlinear eigensolver context
230: .  n     - number of vectors
231: -  is    - set of basis vectors of the initial space

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

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

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

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

246:    Level: intermediate
247: @*/
248: PetscErrorCode NEPSetInitialSpace(NEP nep,PetscInt n,Vec *is)
249: {

255:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)nep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
256:   SlepcBasisReference_Private(n,is,&nep->nini,&nep->IS);
257:   if (n>0) nep->state = NEP_STATE_INITIAL;
258:   return(0);
259: }

263: /*
264:   NEPSetDimensions_Default - Set reasonable values for ncv, mpd if not set
265:   by the user. This is called at setup.
266:  */
267: PetscErrorCode NEPSetDimensions_Default(NEP nep,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
268: {
270:   if (*ncv) { /* ncv set */
271:     if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)nep),1,"The value of ncv must be at least nev");
272:   } else if (*mpd) { /* mpd set */
273:     *ncv = PetscMin(nep->n,nev+(*mpd));
274:   } else { /* neither set: defaults depend on nev being small or large */
275:     if (nev<500) *ncv = PetscMin(nep->n,PetscMax(2*nev,nev+15));
276:     else {
277:       *mpd = 500;
278:       *ncv = PetscMin(nep->n,nev+(*mpd));
279:     }
280:   }
281:   if (!*mpd) *mpd = *ncv;
282:   return(0);
283: }

287: /*@
288:    NEPAllocateSolution - Allocate memory storage for common variables such
289:    as eigenvalues and eigenvectors.

291:    Collective on NEP

293:    Input Parameters:
294: +  nep   - eigensolver context
295: -  extra - number of additional positions, used for methods that require a
296:            working basis slightly larger than ncv

298:    Developers Note:
299:    This is PETSC_EXTERN because it may be required by user plugin NEP
300:    implementations.

302:    Level: developer
303: @*/
304: PetscErrorCode NEPAllocateSolution(NEP nep,PetscInt extra)
305: {
307:   PetscInt       oldsize,newc,requested;
308:   PetscLogDouble cnt;
309:   Mat            T;
310:   Vec            t;

313:   requested = nep->ncv + extra;

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

319:   /* allocate space for eigenvalues and friends */
320:   if (requested != oldsize || !nep->eigr) {
321:     if (oldsize) {
322:       PetscFree4(nep->eigr,nep->eigi,nep->errest,nep->perm);
323:     }
324:     PetscMalloc4(requested,&nep->eigr,requested,&nep->eigi,requested,&nep->errest,requested,&nep->perm);
325:     cnt = newc*sizeof(PetscScalar) + newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
326:     PetscLogObjectMemory((PetscObject)nep,cnt);
327:   }

329:   /* allocate V */
330:   if (!nep->V) { NEPGetBV(nep,&nep->V); }
331:   if (!oldsize) {
332:     if (!((PetscObject)(nep->V))->type_name) {
333:       BVSetType(nep->V,BVSVEC);
334:     }
335:     if (nep->fui==NEP_USER_INTERFACE_SPLIT) T = nep->A[0];
336:     else {
337:       NEPGetFunction(nep,&T,NULL,NULL,NULL);
338:     }
339:     MatCreateVecs(T,&t,NULL);
340:     BVSetSizesFromVec(nep->V,t,requested);
341:     VecDestroy(&t);
342:   } else {
343:     BVResize(nep->V,requested,PETSC_FALSE);
344:   }
345:   return(0);
346: }