Actual source code: qepsetup.c

  1: /*
  2:       QEP routines related to problem setup.

  4:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  5:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  6:    Copyright (c) 2002-2013, 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/qepimpl.h>       /*I "slepcqep.h" I*/
 25: #include <slepc-private/ipimpl.h>

 29: /*@
 30:    QEPSetUp - Sets up all the internal data structures necessary for the
 31:    execution of the QEP solver.

 33:    Collective on QEP

 35:    Input Parameter:
 36: .  qep   - solver context

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

 43:    Level: advanced

 45: .seealso: QEPCreate(), QEPSolve(), QEPDestroy()
 46: @*/
 47: PetscErrorCode QEPSetUp(QEP qep)
 48: {
 50:   PetscBool      khas,mhas,islinear,flg;
 51:   PetscReal      knorm,mnorm;
 52:   Mat            mat[3];

 56:   if (qep->setupcalled) return(0);
 57:   PetscLogEventBegin(QEP_SetUp,qep,0,0,0);

 59:   /* reset the convergence flag from the previous solves */
 60:   qep->reason = QEP_CONVERGED_ITERATING;

 62:   /* Set default solver type (QEPSetFromOptions was not called) */
 63:   if (!((PetscObject)qep)->type_name) {
 64:     QEPSetType(qep,QEPLINEAR);
 65:   }
 66:   PetscObjectTypeCompare((PetscObject)qep,QEPLINEAR,&islinear);
 67:   if (!islinear) {
 68:     if (!qep->st) { QEPGetST(qep,&qep->st); }
 69:     if (!((PetscObject)qep->st)->type_name) {
 70:       STSetType(qep->st,STSHIFT);
 71:     }
 72:   }
 73:   if (!qep->ip) { QEPGetIP(qep,&qep->ip); }
 74:   if (!((PetscObject)qep->ip)->type_name) {
 75:     IPSetType_Default(qep->ip);
 76:   }
 77:   if (!qep->ds) { QEPGetDS(qep,&qep->ds); }
 78:   DSReset(qep->ds);
 79:   if (!((PetscObject)qep->rand)->type_name) {
 80:     PetscRandomSetFromOptions(qep->rand);
 81:   }

 83:   /* Check matrices, transfer them to ST */
 84:   if (!qep->M || !qep->C || !qep->K) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_WRONGSTATE,"QEPSetOperators must be called first");
 85:   if (!islinear) {
 86:     mat[0] = qep->K;
 87:     mat[1] = qep->C;
 88:     mat[2] = qep->M;
 89:     STSetOperators(qep->st,3,mat);
 90:   }

 92:   /* Set problem dimensions */
 93:   MatGetSize(qep->M,&qep->n,NULL);
 94:   MatGetLocalSize(qep->M,&qep->nloc,NULL);
 95:   VecDestroy(&qep->t);
 96:   SlepcMatGetVecsTemplate(qep->M,&qep->t,NULL);
 97:   PetscLogObjectParent(qep,qep->t);

 99:   /* Set default problem type */
100:   if (!qep->problem_type) {
101:     QEPSetProblemType(qep,QEP_GENERAL);
102:   }

104:   /* Compute scaling factor if not set by user */
105:   if (qep->sfactor==0.0) {
106:     MatHasOperation(qep->K,MATOP_NORM,&khas);
107:     MatHasOperation(qep->M,MATOP_NORM,&mhas);
108:     if (khas && mhas) {
109:       MatNorm(qep->K,NORM_INFINITY,&knorm);
110:       MatNorm(qep->M,NORM_INFINITY,&mnorm);
111:       qep->sfactor = PetscSqrtReal(knorm/mnorm);
112:     } else qep->sfactor = 1.0;
113:   }

115:   /* Call specific solver setup */
116:   (*qep->ops->setup)(qep);

118:   /* set tolerance if not yet set */
119:   if (qep->tol==PETSC_DEFAULT) qep->tol = SLEPC_DEFAULT_TOL;

121:   /* set eigenvalue comparison */
122:   switch (qep->which) {
123:     case QEP_LARGEST_MAGNITUDE:
124:       qep->comparison    = SlepcCompareLargestMagnitude;
125:       qep->comparisonctx = NULL;
126:       break;
127:     case QEP_SMALLEST_MAGNITUDE:
128:       qep->comparison    = SlepcCompareSmallestMagnitude;
129:       qep->comparisonctx = NULL;
130:       break;
131:     case QEP_LARGEST_REAL:
132:       qep->comparison    = SlepcCompareLargestReal;
133:       qep->comparisonctx = NULL;
134:       break;
135:     case QEP_SMALLEST_REAL:
136:       qep->comparison    = SlepcCompareSmallestReal;
137:       qep->comparisonctx = NULL;
138:       break;
139:     case QEP_LARGEST_IMAGINARY:
140:       qep->comparison    = SlepcCompareLargestImaginary;
141:       qep->comparisonctx = NULL;
142:       break;
143:     case QEP_SMALLEST_IMAGINARY:
144:       qep->comparison    = SlepcCompareSmallestImaginary;
145:       qep->comparisonctx = NULL;
146:       break;
147:     case QEP_TARGET_MAGNITUDE:
148:       qep->comparison    = SlepcCompareTargetMagnitude;
149:       qep->comparisonctx = &qep->target;
150:       break;
151:     case QEP_TARGET_REAL:
152:       qep->comparison    = SlepcCompareTargetReal;
153:       qep->comparisonctx = &qep->target;
154:       break;
155:     case QEP_TARGET_IMAGINARY:
156:       qep->comparison    = SlepcCompareTargetImaginary;
157:       qep->comparisonctx = &qep->target;
158:       break;
159:   }

161:   if (qep->ncv > 2*qep->n) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_OUTOFRANGE,"ncv must be twice the problem size at most");
162:   if (qep->nev > qep->ncv) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_OUTOFRANGE,"nev bigger than ncv");

164:   /* Setup ST */
165:   if (!islinear) {
166:     PetscObjectTypeCompareAny((PetscObject)qep->st,&flg,STSHIFT,STSINVERT,"");
167:     if (!flg) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_SUP,"Only STSHIFT and STSINVERT spectral transformations can be used in QEP");
168:     STSetUp(qep->st);
169:   }

171:   /* process initial vectors */
172:   if (qep->nini<0) {
173:     qep->nini = -qep->nini;
174:     if (qep->nini>qep->ncv) SETERRQ(PetscObjectComm((PetscObject)qep),1,"The number of initial vectors is larger than ncv");
175:     IPOrthonormalizeBasis_Private(qep->ip,&qep->nini,&qep->IS,qep->V);
176:   }
177:   if (qep->ninil<0) {
178:     if (!qep->leftvecs) { PetscInfo(qep,"Ignoring initial left vectors\n"); }
179:     else {
180:       qep->ninil = -qep->ninil;
181:       if (qep->ninil>qep->ncv) SETERRQ(PetscObjectComm((PetscObject)qep),1,"The number of initial left vectors is larger than ncv");
182:       IPOrthonormalizeBasis_Private(qep->ip,&qep->ninil,&qep->ISL,qep->W);
183:     }
184:   }
185:   PetscLogEventEnd(QEP_SetUp,qep,0,0,0);
186:   qep->setupcalled = 1;
187:   return(0);
188: }

192: /*@
193:    QEPSetOperators - Sets the matrices associated with the quadratic eigenvalue problem.

195:    Collective on QEP and Mat

197:    Input Parameters:
198: +  qep - the eigenproblem solver context
199: .  M   - the first coefficient matrix
200: .  C   - the second coefficient matrix
201: -  K   - the third coefficient matrix

203:    Notes:
204:    The quadratic eigenproblem is defined as (l^2*M + l*C + K)*x = 0, where l is
205:    the eigenvalue and x is the eigenvector.

207:    Level: beginner

209: .seealso: QEPSolve(), QEPGetOperators()
210: @*/
211: PetscErrorCode QEPSetOperators(QEP qep,Mat M,Mat C,Mat K)
212: {
214:   PetscInt       m,n,m0;


225:   /* Check for square matrices */
226:   MatGetSize(M,&m,&n);
227:   if (m!=n) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_WRONG,"M is a non-square matrix");
228:   m0=m;
229:   MatGetSize(C,&m,&n);
230:   if (m!=n) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_WRONG,"C is a non-square matrix");
231:   if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_INCOMP,"Dimensions of M and C do not match");
232:   MatGetSize(K,&m,&n);
233:   if (m!=n) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_WRONG,"K is a non-square matrix");
234:   if (m!=m0) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_INCOMP,"Dimensions of M and K do not match");

236:   /* Store a copy of the matrices */
237:   if (qep->setupcalled) { QEPReset(qep); }
238:   PetscObjectReference((PetscObject)M);
239:   MatDestroy(&qep->M);
240:   qep->M = M;
241:   PetscObjectReference((PetscObject)C);
242:   MatDestroy(&qep->C);
243:   qep->C = C;
244:   PetscObjectReference((PetscObject)K);
245:   MatDestroy(&qep->K);
246:   qep->K = K;
247:   return(0);
248: }

252: /*@
253:    QEPGetOperators - Gets the matrices associated with the quadratic eigensystem.

255:    Collective on QEP and Mat

257:    Input Parameter:
258: .  qep - the QEP context

260:    Output Parameters:
261: +  M   - the first coefficient matrix
262: .  C   - the second coefficient matrix
263: -  K   - the third coefficient matrix

265:    Level: intermediate

267: .seealso: QEPSolve(), QEPSetOperators()
268: @*/
269: PetscErrorCode QEPGetOperators(QEP qep,Mat *M,Mat *C,Mat *K)
270: {
276:   return(0);
277: }

281: /*@
282:    QEPSetInitialSpace - Specify a basis of vectors that constitute the initial
283:    space, that is, the subspace from which the solver starts to iterate.

285:    Collective on QEP and Vec

287:    Input Parameter:
288: +  qep   - the quadratic eigensolver context
289: .  n     - number of vectors
290: -  is    - set of basis vectors of the initial space

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

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

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

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

305:    Level: intermediate

307: .seealso: QEPSetInitialSpaceLeft()
308: @*/
309: PetscErrorCode QEPSetInitialSpace(QEP qep,PetscInt n,Vec *is)
310: {

316:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
317:   SlepcBasisReference_Private(n,is,&qep->nini,&qep->IS);
318:   if (n>0) qep->setupcalled = 0;
319:   return(0);
320: }

324: /*@
325:    QEPSetInitialSpaceLeft - Specify a basis of vectors that constitute the initial
326:    left space, that is, the subspace from which the solver starts to iterate for
327:    building the left subspace (in methods that work with two subspaces).

329:    Collective on QEP and Vec

331:    Input Parameter:
332: +  qep   - the quadratic eigensolver context
333: .  n     - number of vectors
334: -  is    - set of basis vectors of the initial left space

336:    Notes:
337:    Some solvers start to iterate on a single vector (initial left vector). In that case,
338:    the other vectors are ignored.

340:    These vectors do not persist from one QEPSolve() call to the other, so the
341:    initial left space should be set every time.

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

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

349:    Level: intermediate

351: .seealso: QEPSetInitialSpace()
352: @*/
353: PetscErrorCode QEPSetInitialSpaceLeft(QEP qep,PetscInt n,Vec *is)
354: {

360:   if (n<0) SETERRQ(PetscObjectComm((PetscObject)qep),PETSC_ERR_ARG_OUTOFRANGE,"Argument n cannot be negative");
361:   SlepcBasisReference_Private(n,is,&qep->ninil,&qep->ISL);
362:   if (n>0) qep->setupcalled = 0;
363:   return(0);
364: }

368: /*
369:   QEPAllocateSolution - Allocate memory storage for common variables such
370:   as eigenvalues and eigenvectors. All vectors in V (and W) share a
371:   contiguous chunk of memory.
372: */
373: PetscErrorCode QEPAllocateSolution(QEP qep)
374: {
376:   PetscInt       newc,cnt;

379:   if (qep->allocated_ncv != qep->ncv) {
380:     newc = PetscMax(0,qep->ncv-qep->allocated_ncv);
381:     QEPFreeSolution(qep);
382:     cnt = 0;
383:     PetscMalloc(qep->ncv*sizeof(PetscScalar),&qep->eigr);
384:     PetscMalloc(qep->ncv*sizeof(PetscScalar),&qep->eigi);
385:     cnt += 2*newc*sizeof(PetscScalar);
386:     PetscMalloc(qep->ncv*sizeof(PetscReal),&qep->errest);
387:     PetscMalloc(qep->ncv*sizeof(PetscInt),&qep->perm);
388:     cnt += 2*newc*sizeof(PetscReal);
389:     PetscLogObjectMemory(qep,cnt);
390:     VecDuplicateVecs(qep->t,qep->ncv,&qep->V);
391:     PetscLogObjectParents(qep,qep->ncv,qep->V);
392:     qep->allocated_ncv = qep->ncv;
393:   }
394:   return(0);
395: }

399: /*
400:   QEPFreeSolution - Free memory storage. This routine is related to
401:   QEPAllocateSolution().
402: */
403: PetscErrorCode QEPFreeSolution(QEP qep)
404: {

408:   if (qep->allocated_ncv > 0) {
409:     PetscFree(qep->eigr);
410:     PetscFree(qep->eigi);
411:     PetscFree(qep->errest);
412:     PetscFree(qep->perm);
413:     VecDestroyVecs(qep->allocated_ncv,&qep->V);
414:     qep->allocated_ncv = 0;
415:   }
416:   return(0);
417: }