Actual source code: setup.c

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

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

  8:    This file is part of SLEPc.
  9:       
 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 private/epsimpl.h

 28: /*@
 29:    EPSSetUp - Sets up all the internal data structures necessary for the
 30:    execution of the eigensolver. Then calls STSetUp() for any set-up
 31:    operations associated to the ST object.

 33:    Collective on EPS

 35:    Input Parameter:
 36: .  eps   - eigenproblem solver context

 38:    Level: advanced

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

 45:    This function sets a random initial vector if none has been provided.

 47: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp()
 48: @*/
 49: PetscErrorCode EPSSetUp(EPS eps)
 50: {
 52:   PetscInt       i;
 53:   Vec            v0,w0;
 54:   Mat            A,B;
 55:   PetscInt       N;
 56: 

 60:   if (eps->setupcalled) return(0);

 62:   PetscLogEventBegin(EPS_SetUp,eps,0,0,0);

 64:   /* Set default solver type */
 65:   if (!((PetscObject)eps)->type_name) {
 66:     EPSSetType(eps,EPSKRYLOVSCHUR);
 67:   }
 68: 
 69:   STGetOperators(eps->OP,&A,&B);
 70:   /* Set default problem type */
 71:   if (!eps->problem_type) {
 72:     if (B==PETSC_NULL) {
 73:       EPSSetProblemType(eps,EPS_NHEP);
 74:     }
 75:     else {
 76:       EPSSetProblemType(eps,EPS_GNHEP);
 77:     }
 78:   } else if ((B && !eps->isgeneralized) || (!B && eps->isgeneralized)) {
 79:     SETERRQ(0,"Warning: Inconsistent EPS state");
 80:   }
 81: 
 82:   if (eps->ispositive) {
 83:     STGetBilinearForm(eps->OP,&B);
 84:     IPSetBilinearForm(eps->ip,B,IPINNER_HERMITIAN);
 85:     MatDestroy(B);
 86:   } else {
 87:     IPSetBilinearForm(eps->ip,PETSC_NULL,IPINNER_HERMITIAN);
 88:   }
 89: 
 90:   /* Create random initial vectors if not set */
 91:   /* right */
 92:   EPSGetInitialVector(eps,&v0);
 93:   if (!v0) {
 94:     MatGetVecs(A,&v0,PETSC_NULL);
 95:     SlepcVecSetRandom(v0);
 96:     eps->vec_initial = v0;
 97:   }
 98:   /* left */
 99:   EPSGetLeftInitialVector(eps,&w0);
100:   if (!w0) {
101:     MatGetVecs(A,PETSC_NULL,&w0);
102:     SlepcVecSetRandom(w0);
103:     eps->vec_initial_left = w0;
104:   }

106:   VecGetSize(eps->vec_initial,&N);
107:   if (eps->nev > N) eps->nev = N;
108:   if (eps->ncv > N) eps->ncv = N;

110:   (*eps->ops->setup)(eps);
111:   STSetUp(eps->OP);

113:   /* DSV is equal to the columns of DS followed by the ones in V */
114:   PetscFree(eps->DSV);
115:   PetscMalloc((eps->ncv+eps->nds)*sizeof(Vec),&eps->DSV);
116:   for (i = 0; i < eps->nds; i++) eps->DSV[i] = eps->DS[i];
117:   for (i = 0; i < eps->ncv; i++) eps->DSV[i+eps->nds] = eps->V[i];
118: 
119:   if (eps->nds>0) {
120:     if (!eps->ds_ortho) {
121:       /* orthonormalize vectors in DS if necessary */
122:       IPQRDecomposition(eps->ip,eps->DS,0,eps->nds,PETSC_NULL,0,PETSC_NULL);
123:     }
124:     IPOrthogonalize(eps->ip,eps->nds,PETSC_NULL,eps->DS,eps->vec_initial,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL,PETSC_NULL);
125:   }

127:   STCheckNullSpace(eps->OP,eps->nds,eps->DS);
128: 
129:   PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
130:   eps->setupcalled = 1;
131:   return(0);
132: }

136: /*@
137:    EPSSetInitialVector - Sets the initial vector from which the 
138:    eigensolver starts to iterate.

140:    Collective on EPS and Vec

142:    Input Parameters:
143: +  eps - the eigensolver context
144: -  vec - the vector

146:    Level: intermediate

148: .seealso: EPSGetInitialVector(), EPSSetLeftInitialVector()

150: @*/
151: PetscErrorCode EPSSetInitialVector(EPS eps,Vec vec)
152: {
154: 
159:   PetscObjectReference((PetscObject)vec);
160:   if (eps->vec_initial) {
161:     VecDestroy(eps->vec_initial);
162:   }
163:   eps->vec_initial = vec;
164:   return(0);
165: }

169: /*@
170:    EPSGetInitialVector - Gets the initial vector associated with the 
171:    eigensolver; if the vector was not set it will return a 0 pointer or
172:    a vector randomly generated by EPSSetUp().

174:    Not collective, but vector is shared by all processors that share the EPS

176:    Input Parameter:
177: .  eps - the eigensolver context

179:    Output Parameter:
180: .  vec - the vector

182:    Level: intermediate

184: .seealso: EPSSetInitialVector(), EPSGetLeftInitialVector()

186: @*/
187: PetscErrorCode EPSGetInitialVector(EPS eps,Vec *vec)
188: {
192:   *vec = eps->vec_initial;
193:   return(0);
194: }

198: /*@
199:    EPSSetLeftInitialVector - Sets the initial vector from which the eigensolver 
200:    starts to iterate, corresponding to the left recurrence (two-sided solvers).

202:    Collective on EPS and Vec

204:    Input Parameters:
205: +  eps - the eigensolver context
206: -  vec - the vector

208:    Level: intermediate

210: .seealso: EPSGetLeftInitialVector(), EPSSetInitialVector()

212: @*/
213: PetscErrorCode EPSSetLeftInitialVector(EPS eps,Vec vec)
214: {
216: 
221:   PetscObjectReference((PetscObject)vec);
222:   if (eps->vec_initial_left) {
223:     VecDestroy(eps->vec_initial_left);
224:   }
225:   eps->vec_initial_left = vec;
226:   return(0);
227: }

231: /*@
232:    EPSGetLeftInitialVector - Gets the left initial vector associated with the 
233:    eigensolver; if the vector was not set it will return a 0 pointer or
234:    a vector randomly generated by EPSSetUp().

236:    Not collective, but vector is shared by all processors that share the EPS

238:    Input Parameter:
239: .  eps - the eigensolver context

241:    Output Parameter:
242: .  vec - the vector

244:    Level: intermediate

246: .seealso: EPSSetLeftInitialVector(), EPSGetLeftInitialVector()

248: @*/
249: PetscErrorCode EPSGetLeftInitialVector(EPS eps,Vec *vec)
250: {
254:   *vec = eps->vec_initial_left;
255:   return(0);
256: }

260: /*@
261:    EPSSetOperators - Sets the matrices associated with the eigenvalue problem.

263:    Collective on EPS and Mat

265:    Input Parameters:
266: +  eps - the eigenproblem solver context
267: .  A  - the matrix associated with the eigensystem
268: -  B  - the second matrix in the case of generalized eigenproblems

270:    Notes: 
271:    To specify a standard eigenproblem, use PETSC_NULL for parameter B.

273:    Level: beginner

275: .seealso: EPSSolve(), EPSGetST(), STGetOperators()
276: @*/
277: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
278: {
280:   PetscInt       m,n;


289:   /* Check for square matrices */
290:   MatGetSize(A,&m,&n);
291:   if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
292:   if (B) {
293:     MatGetSize(B,&m,&n);
294:     if (m!=n) { SETERRQ(1,"B is a non-square matrix"); }
295:   }

297:   STSetOperators(eps->OP,A,B);
298:   eps->setupcalled = 0;  /* so that next solve call will call setup */

300:   /* Destroy randomly generated initial vectors */
301:   if (eps->vec_initial) {
302:     VecDestroy(eps->vec_initial);
303:     eps->vec_initial = PETSC_NULL;
304:   }
305:   if (eps->vec_initial_left) {
306:     VecDestroy(eps->vec_initial_left);
307:     eps->vec_initial_left = PETSC_NULL;
308:   }

310:   return(0);
311: }

315: /*@
316:    EPSGetOperators - Gets the matrices associated with the eigensystem.

318:    Collective on EPS and Mat

320:    Input Parameter:
321: .  eps - the EPS context

323:    Output Parameters:
324: +  A  - the matrix associated with the eigensystem
325: -  B  - the second matrix in the case of generalized eigenproblems

327:    Level: intermediate

329: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
330: @*/
331: PetscErrorCode EPSGetOperators(EPS eps, Mat *A, Mat *B)
332: {
334:   ST             st;

340:   EPSGetST(eps,&st);
341:   STGetOperators(st,A,B);
342:   return(0);
343: }

347: /*@
348:    EPSAttachDeflationSpace - Add vectors to the basis of the deflation space.

350:    Not Collective

352:    Input Parameter:
353: +  eps   - the eigenproblem solver context
354: .  n     - number of vectors to add
355: .  ds    - set of basis vectors of the deflation space
356: -  ortho - PETSC_TRUE if basis vectors of deflation space are orthonormal

358:    Notes:
359:    When a deflation space is given, the eigensolver seeks the eigensolution
360:    in the restriction of the problem to the orthogonal complement of this
361:    space. This can be used for instance in the case that an invariant 
362:    subspace is known beforehand (such as the nullspace of the matrix).

364:    The basis vectors can be provided all at once or incrementally with
365:    several calls to EPSAttachDeflationSpace().

367:    Use a value of PETSC_TRUE for parameter ortho if all the vectors passed
368:    in are known to be mutually orthonormal.

370:    Level: intermediate

372: .seealso: EPSRemoveDeflationSpace()
373: @*/
374: PetscErrorCode EPSAttachDeflationSpace(EPS eps,PetscInt n,Vec *ds,PetscTruth ortho)
375: {
377:   PetscInt       i;
378:   Vec            *tvec;
379: 
382:   tvec = eps->DS;
383:   if (n+eps->nds > 0) {
384:      PetscMalloc((n+eps->nds)*sizeof(Vec), &eps->DS);
385:   }
386:   if (eps->nds > 0) {
387:     for (i=0; i<eps->nds; i++) eps->DS[i] = tvec[i];
388:     PetscFree(tvec);
389:   }
390:   for (i=0; i<n; i++) {
391:     VecDuplicate(ds[i],&eps->DS[i + eps->nds]);
392:     VecCopy(ds[i],eps->DS[i + eps->nds]);
393:   }
394:   eps->nds += n;
395:   if (!ortho) eps->ds_ortho = PETSC_FALSE;
396:   eps->setupcalled = 0;
397:   return(0);
398: }

402: /*@
403:    EPSRemoveDeflationSpace - Removes the deflation space.

405:    Not Collective

407:    Input Parameter:
408: .  eps   - the eigenproblem solver context

410:    Level: intermediate

412: .seealso: EPSAttachDeflationSpace()
413: @*/
414: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
415: {
417: 
420:   if (eps->nds > 0) {
421:     VecDestroyVecs(eps->DS, eps->nds);
422:   }
423:   eps->ds_ortho = PETSC_TRUE;
424:   eps->setupcalled = 0;
425:   return(0);
426: }