Actual source code: pepsetup.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:       PEP routines related to problem setup.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2014, 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/pepimpl.h>       /*I "slepcpep.h" I*/

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

 32:    Collective on PEP

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

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

 42:    Level: advanced

 44: .seealso: PEPCreate(), PEPSolve(), PEPDestroy()
 45: @*/
 46: PetscErrorCode PEPSetUp(PEP pep)
 47: {
 49:   SlepcSC        sc;
 50:   PetscBool      islinear,istrivial,flg;
 51:   PetscInt       i,k;

 55:   if (pep->state) return(0);
 56:   PetscLogEventBegin(PEP_SetUp,pep,0,0,0);

 58:   /* reset the convergence flag from the previous solves */
 59:   pep->reason = PEP_CONVERGED_ITERATING;

 61:   /* set default solver type (PEPSetFromOptions was not called) */
 62:   if (!((PetscObject)pep)->type_name) {
 63:     PEPSetType(pep,PEPTOAR);
 64:   }
 65:   if (!pep->st) { PEPGetST(pep,&pep->st); }
 66:   PetscObjectTypeCompare((PetscObject)pep,PEPLINEAR,&islinear);
 67:   if (!islinear) {
 68:     if (!((PetscObject)pep->st)->type_name) {
 69:       STSetType(pep->st,STSHIFT);
 70:     }
 71:   }
 72:   if (!pep->ds) { PEPGetDS(pep,&pep->ds); }
 73:   DSReset(pep->ds);
 74:   if (!pep->rg) { PEPGetRG(pep,&pep->rg); }
 75:   if (!((PetscObject)pep->rg)->type_name) {
 76:     RGSetType(pep->rg,RGINTERVAL);
 77:   }
 78:   if (!((PetscObject)pep->rand)->type_name) {
 79:     PetscRandomSetFromOptions(pep->rand);
 80:   }

 82:   /* check matrices, transfer them to ST */
 83:   if (!pep->A) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONGSTATE,"PEPSetOperators must be called first");
 84:   STSetOperators(pep->st,pep->nmat,pep->A);

 86:   /* set problem dimensions */
 87:   MatGetSize(pep->A[0],&pep->n,NULL);
 88:   MatGetLocalSize(pep->A[0],&pep->nloc,NULL);

 90:   /* set default problem type */
 91:   if (!pep->problem_type) {
 92:     PEPSetProblemType(pep,PEP_GENERAL);
 93:   }

 95:   /* initialization of matrix norms */
 96:   if (pep->conv==PEP_CONV_NORM) {
 97:     for (i=0;i<pep->nmat;i++) {
 98:       if (!pep->nrma[i]) {
 99:         MatNorm(pep->A[i],NORM_INFINITY,&pep->nrma[i]);
100:       }
101:     }
102:   }

104:   /* call specific solver setup */
105:   (*pep->ops->setup)(pep);

107:   /* set tolerance if not yet set */
108:   if (pep->tol==PETSC_DEFAULT) pep->tol = SLEPC_DEFAULT_TOL;
109:   if (pep->refine) {
110:     if (pep->rtol==PETSC_DEFAULT) pep->rtol = SLEPC_DEFAULT_TOL;
111:     if (pep->rits==PETSC_DEFAULT) pep->rits = (pep->refine==PEP_REFINE_SIMPLE)? 10: 1;
112:   }

114:   /* fill sorting criterion context */
115:   switch (pep->which) {
116:     case PEP_LARGEST_MAGNITUDE:
117:       pep->sc->comparison    = SlepcCompareLargestMagnitude;
118:       pep->sc->comparisonctx = NULL;
119:       break;
120:     case PEP_SMALLEST_MAGNITUDE:
121:       pep->sc->comparison    = SlepcCompareSmallestMagnitude;
122:       pep->sc->comparisonctx = NULL;
123:       break;
124:     case PEP_LARGEST_REAL:
125:       pep->sc->comparison    = SlepcCompareLargestReal;
126:       pep->sc->comparisonctx = NULL;
127:       break;
128:     case PEP_SMALLEST_REAL:
129:       pep->sc->comparison    = SlepcCompareSmallestReal;
130:       pep->sc->comparisonctx = NULL;
131:       break;
132:     case PEP_LARGEST_IMAGINARY:
133:       pep->sc->comparison    = SlepcCompareLargestImaginary;
134:       pep->sc->comparisonctx = NULL;
135:       break;
136:     case PEP_SMALLEST_IMAGINARY:
137:       pep->sc->comparison    = SlepcCompareSmallestImaginary;
138:       pep->sc->comparisonctx = NULL;
139:       break;
140:     case PEP_TARGET_MAGNITUDE:
141:       pep->sc->comparison    = SlepcCompareTargetMagnitude;
142:       pep->sc->comparisonctx = &pep->target;
143:       break;
144:     case PEP_TARGET_REAL:
145:       pep->sc->comparison    = SlepcCompareTargetReal;
146:       pep->sc->comparisonctx = &pep->target;
147:       break;
148:     case PEP_TARGET_IMAGINARY:
149:       pep->sc->comparison    = SlepcCompareTargetImaginary;
150:       pep->sc->comparisonctx = &pep->target;
151:       break;
152:   }
153:   pep->sc->map    = NULL;
154:   pep->sc->mapobj = NULL;

156:   /* fill sorting criterion for DS */
157:   DSGetSlepcSC(pep->ds,&sc);
158:   RGIsTrivial(pep->rg,&istrivial);
159:   sc->rg            = istrivial? NULL: pep->rg;
160:   sc->comparison    = pep->sc->comparison;
161:   sc->comparisonctx = pep->sc->comparisonctx;
162:   sc->map           = SlepcMap_ST;
163:   sc->mapobj        = (PetscObject)pep->st;

165:   /* setup ST */
166:   if (!islinear) {
167:     PetscObjectTypeCompareAny((PetscObject)pep->st,&flg,STSHIFT,STSINVERT,"");
168:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used in PEP");
169:     STSetUp(pep->st);
170:     /* compute matrix coefficients */
171:     STGetTransform(pep->st,&flg);
172:     if (!flg) {
173:       STMatSetUp(pep->st,1.0,pep->solvematcoeffs);
174:     } else {
175:       if (pep->basis!=PEP_BASIS_MONOMIAL) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_SUP,"Cannot use ST-transform with non-monomial basis in PEP");
176:     }
177:   }

179:   /* compute scale factor if no set by user */
180:   PEPComputeScaleFactor(pep);

182:   /* build balancing matrix if required */
183:   if (pep->scale==PEP_SCALE_DIAGONAL || pep->scale==PEP_SCALE_BOTH) {
184:     if (!pep->Dl) {
185:       BVGetVec(pep->V,&pep->Dl);
186:       PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->Dl);
187:     }
188:     if (!pep->Dr) {
189:       BVGetVec(pep->V,&pep->Dr);
190:       PetscLogObjectParent((PetscObject)pep,(PetscObject)pep->Dr);
191:     }
192:     PEPBuildDiagonalScaling(pep);
193:   }

195:   /* process initial vectors */
196:   if (pep->nini<0) {
197:     k = -pep->nini;
198:     if (k>pep->ncv) SETERRQ(PetscObjectComm((PetscObject)pep),1,"The number of initial vectors is larger than ncv");
199:     BVInsertVecs(pep->V,0,&k,pep->IS,PETSC_TRUE);
200:     SlepcBasisDestroy_Private(&pep->nini,&pep->IS);
201:     pep->nini = k;
202:   }
203:   PetscLogEventEnd(PEP_SetUp,pep,0,0,0);
204:   pep->state = PEP_STATE_SETUP;
205:   return(0);
206: }

210: /*@
211:    PEPSetOperators - Sets the coefficient matrices associated with the polynomial
212:    eigenvalue problem.

214:    Collective on PEP and Mat

216:    Input Parameters:
217: +  pep  - the eigenproblem solver context
218: .  nmat - number of matrices in array A
219: -  A    - the array of matrices associated with the eigenproblem

221:    Notes:
222:    The polynomial eigenproblem is defined as P(l)*x=0, where l is
223:    the eigenvalue, x is the eigenvector, and P(l) is defined as
224:    P(l) = A_0 + l*A_1 + ... + l^d*A_d, with d=nmat-1 (the degree of P).
225:    For non-monomial bases, this expression is different.

227:    Level: beginner

229: .seealso: PEPSolve(), PEPGetOperators(), PEPGetNumMatrices(), PEPSetBasis()
230: @*/
231: PetscErrorCode PEPSetOperators(PEP pep,PetscInt nmat,Mat A[])
232: {
234:   PetscInt       i,n,m,m0=0;

239:   if (nmat <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Must have one or more matrices, you have %D",nmat);

242:   if (pep->state) { PEPReset(pep); }
243:   PetscMalloc1(nmat,&pep->A);
244:   PetscCalloc3(3*nmat,&pep->pbc,nmat,&pep->solvematcoeffs,nmat,&pep->nrma);
245:   for (i=0;i<nmat;i++) pep->pbc[i] = 1.0;  /* default to monomial basis */
246:   PetscLogObjectMemory((PetscObject)pep,nmat*sizeof(Mat)+4*nmat*sizeof(PetscReal)+nmat*sizeof(PetscScalar));
247:   for (i=0;i<nmat;i++) {
250:     MatGetSize(A[i],&m,&n);
251:     if (m!=n) SETERRQ1(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_WRONG,"A[%D] is a non-square matrix",i);
252:     if (!i) m0 = m;
253:     if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_INCOMP,"Dimensions of matrices do not match with each other");
254:     PetscObjectReference((PetscObject)A[i]);
255:     pep->A[i] = A[i];
256:   }
257:   pep->nmat = nmat;
258:   return(0);
259: }

263: /*@
264:    PEPGetOperators - Gets the matrices associated with the polynomial eigensystem.

266:    Not collective, though parallel Mats are returned if the PEP is parallel

268:    Input Parameters:
269: +  pep - the PEP context
270: -  k   - the index of the requested matrix (starting in 0)

272:    Output Parameter:
273: .  A - the requested matrix

275:    Level: intermediate

277: .seealso: PEPSolve(), PEPSetOperators(), PEPGetNumMatrices()
278: @*/
279: PetscErrorCode PEPGetOperators(PEP pep,PetscInt k,Mat *A)
280: {
284:   if (k<0 || k>=pep->nmat) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",pep->nmat-1);
285:   *A = pep->A[k];
286:   return(0);
287: }

291: /*@
292:    PEPGetNumMatrices - Returns the number of matrices stored in the PEP.

294:    Not collective

296:    Input Parameter:
297: .  pep - the PEP context

299:    Output Parameters:
300: .  nmat - the number of matrices passed in PEPSetOperators()

302:    Level: intermediate

304: .seealso: PEPSetOperators()
305: @*/
306: PetscErrorCode PEPGetNumMatrices(PEP pep,PetscInt *nmat)
307: {
311:   *nmat = pep->nmat;
312:   return(0);
313: }

317: /*@
318:    PEPSetInitialSpace - Specify a basis of vectors that constitute the initial
319:    space, that is, the subspace from which the solver starts to iterate.

321:    Collective on PEP and Vec

323:    Input Parameter:
324: +  pep   - the polynomial eigensolver context
325: .  n     - number of vectors
326: -  is    - set of basis vectors of the initial space

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

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

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

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

341:    Level: intermediate
342: @*/
343: PetscErrorCode PEPSetInitialSpace(PEP pep,PetscInt n,Vec *is)
344: {

350:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)pep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
351:   SlepcBasisReference_Private(n,is,&pep->nini,&pep->IS);
352:   if (n>0) pep->state = PEP_STATE_INITIAL;
353:   return(0);
354: }

358: /*
359:   PEPSetDimensions_Default - Set reasonable values for ncv, mpd if not set
360:   by the user. This is called at setup.
361:  */
362: PetscErrorCode PEPSetDimensions_Default(PEP pep,PetscInt nev,PetscInt *ncv,PetscInt *mpd)
363: {
365:   PetscBool      krylov;
366:   PetscInt       dim;

369:   PetscObjectTypeCompareAny((PetscObject)pep,&krylov,PEPTOAR,PEPQARNOLDI,"");
370:   dim = krylov?(pep->nmat-1)*pep->n:pep->n;
371:   if (*ncv) { /* ncv set */
372:     if (krylov) {
373:       if (*ncv<nev+1 && !(*ncv==nev && *ncv==dim)) SETERRQ(PetscObjectComm((PetscObject)pep),1,"The value of ncv must be at least nev+1");
374:     } else {
375:       if (*ncv<nev) SETERRQ(PetscObjectComm((PetscObject)pep),1,"The value of ncv must be at least nev");
376:     }
377:   } else if (*mpd) { /* mpd set */
378:     *ncv = PetscMin(dim,nev+(*mpd));
379:   } else { /* neither set: defaults depend on nev being small or large */
380:     if (nev<500) *ncv = PetscMin(dim,PetscMax(2*nev,nev+15));
381:     else {
382:       *mpd = 500;
383:       *ncv = PetscMin(dim,nev+(*mpd));
384:     }
385:   }
386:   if (!*mpd) *mpd = *ncv;
387:   return(0);
388: }

392: /*@
393:    PEPAllocateSolution - Allocate memory storage for common variables such
394:    as eigenvalues and eigenvectors.

396:    Collective on PEP

398:    Input Parameters:
399: +  pep   - eigensolver context
400: -  extra - number of additional positions, used for methods that require a
401:            working basis slightly larger than ncv

403:    Developers Note:
404:    This is PETSC_EXTERN because it may be required by user plugin PEP
405:    implementations.

407:    Level: developer
408: @*/
409: PetscErrorCode PEPAllocateSolution(PEP pep,PetscInt extra)
410: {
412:   PetscInt       oldsize,newc,requested;
413:   PetscLogDouble cnt;
414:   Vec            t;

417:   requested = pep->ncv + extra;

419:   /* oldsize is zero if this is the first time setup is called */
420:   BVGetSizes(pep->V,NULL,NULL,&oldsize);
421:   newc = PetscMax(0,requested-oldsize);

423:   /* allocate space for eigenvalues and friends */
424:   if (requested != oldsize) {
425:     if (oldsize) {
426:       PetscFree4(pep->eigr,pep->eigi,pep->errest,pep->perm);
427:     }
428:     PetscMalloc4(requested,&pep->eigr,requested,&pep->eigi,requested,&pep->errest,requested,&pep->perm);
429:     cnt = 2*newc*sizeof(PetscScalar) + newc*sizeof(PetscReal) + newc*sizeof(PetscInt);
430:     PetscLogObjectMemory((PetscObject)pep,cnt);
431:   }

433:   /* allocate V */
434:   if (!pep->V) { PEPGetBV(pep,&pep->V); }
435:   if (!oldsize) {
436:     if (!((PetscObject)(pep->V))->type_name) {
437:       BVSetType(pep->V,BVSVEC);
438:     }
439:     STMatGetVecs(pep->st,&t,NULL);
440:     BVSetSizesFromVec(pep->V,t,requested);
441:     VecDestroy(&t);
442:   } else {
443:     BVResize(pep->V,requested,PETSC_FALSE);
444:   }
445:   return(0);
446: }