Actual source code: setup.c

  1: /*
  2:       EPS routines related to problem setup.
  3: */
 4:  #include src/eps/epsimpl.h

  8: /*@
  9:    EPSSetUp - Sets up all the internal data structures necessary for the
 10:    execution of the eigensolver. Then calls STSetUp() for any set-up
 11:    operations associated to the ST object.

 13:    Collective on EPS

 15:    Input Parameter:
 16: .  eps   - eigenproblem solver context

 18:    Level: advanced

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

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

 27: .seealso: EPSCreate(), EPSSolve(), EPSDestroy(), STSetUp()
 28: @*/
 29: PetscErrorCode EPSSetUp(EPS eps)
 30: {
 32:   int            i;
 33:   Vec            v0,w0;
 34:   Mat            A,B;
 35: 

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

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

 43:   /* Set default solver type */
 44:   if (!eps->type_name) {
 45:     EPSSetType(eps,EPSARNOLDI);
 46:   }

 48:   /* Set default eta for orthogonalization */
 49:   if (eps->orthog_eta == PETSC_DEFAULT)
 50:     eps->orthog_eta = 0.7071;
 51: 
 52:   STGetOperators(eps->OP,&A,&B);
 53:   /* Set default problem type */
 54:   if (!eps->problem_type) {
 55:     if (B==PETSC_NULL) {
 56:       EPSSetProblemType(eps,EPS_NHEP);
 57:     }
 58:     else {
 59:       EPSSetProblemType(eps,EPS_GNHEP);
 60:     }
 61:   } else if ((B && !eps->isgeneralized) || (!B && eps->isgeneralized)) {
 62:     SETERRQ(0,"Warning: Inconsistent EPS state");
 63:   }
 64: 
 65:   /* Create random initial vectors if not set */
 66:   /* right */
 67:   EPSGetInitialVector(eps,&v0);
 68:   if (!v0) {
 69:     MatGetVecs(A,&v0,PETSC_NULL);
 70:     SlepcVecSetRandom(v0);
 71:     eps->vec_initial = v0;
 72:   }
 73:   /* left */
 74:   EPSGetLeftInitialVector(eps,&w0);
 75:   if (!w0) {
 76:     MatGetVecs(A,PETSC_NULL,&w0);
 77:     SlepcVecSetRandom(w0);
 78:     eps->vec_initial_left = w0;
 79:   }

 81:   (*eps->ops->setup)(eps);
 82:   STSetUp(eps->OP);

 84:   /* DSV is equal to the columns of DS followed by the ones in V */
 85:   PetscFree(eps->DSV);
 86:   PetscMalloc((eps->ncv+eps->nds)*sizeof(Vec),&eps->DSV);
 87:   for (i = 0; i < eps->nds; i++) eps->DSV[i] = eps->DS[i];
 88:   for (i = 0; i < eps->ncv; i++) eps->DSV[i+eps->nds] = eps->V[i];
 89: 
 90:   if (eps->nds>0) {
 91:     if (!eps->ds_ortho) {
 92:       /* orthonormalize vectors in DS if necessary */
 93:       EPSQRDecomposition(eps,eps->DS,0,eps->nds,PETSC_NULL,0);
 94:     }
 95:     EPSOrthogonalize(eps,eps->nds,eps->DS,eps->vec_initial,PETSC_NULL,PETSC_NULL,PETSC_NULL);
 96:   }

 98:   STCheckNullSpace(eps->OP,eps->nds,eps->DS);
 99: 
100:   PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
101:   eps->setupcalled = 1;
102:   return(0);
103: }

107: /*@
108:    EPSSetInitialVector - Sets the initial vector from which the 
109:    eigensolver starts to iterate.

111:    Collective on EPS and Vec

113:    Input Parameters:
114: +  eps - the eigensolver context
115: -  vec - the vector

117:    Level: intermediate

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

121: @*/
122: PetscErrorCode EPSSetInitialVector(EPS eps,Vec vec)
123: {
125: 
130:   if (eps->vec_initial) {
131:     VecDestroy(eps->vec_initial);
132:   }
133:   eps->vec_initial = vec;
134:   PetscObjectReference((PetscObject)eps->vec_initial);
135:   return(0);
136: }

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

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

147:    Input Parameter:
148: .  eps - the eigensolver context

150:    Output Parameter:
151: .  vec - the vector

153:    Level: intermediate

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

157: @*/
158: PetscErrorCode EPSGetInitialVector(EPS eps,Vec *vec)
159: {
162:   *vec = eps->vec_initial;
163:   return(0);
164: }

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

172:    Collective on EPS and Vec

174:    Input Parameters:
175: +  eps - the eigensolver context
176: -  vec - the vector

178:    Level: intermediate

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

182: @*/
183: PetscErrorCode EPSSetLeftInitialVector(EPS eps,Vec vec)
184: {
186: 
191:   if (eps->vec_initial_left) {
192:     VecDestroy(eps->vec_initial_left);
193:   }
194:   eps->vec_initial_left = vec;
195:   PetscObjectReference((PetscObject)eps->vec_initial_left);
196:   return(0);
197: }

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

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

208:    Input Parameter:
209: .  eps - the eigensolver context

211:    Output Parameter:
212: .  vec - the vector

214:    Level: intermediate

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

218: @*/
219: PetscErrorCode EPSGetLeftInitialVector(EPS eps,Vec *vec)
220: {
223:   *vec = eps->vec_initial_left;
224:   return(0);
225: }

229: /*@
230:    EPSSetOperators - Sets the matrices associated with the eigenvalue problem.

232:    Collective on EPS and Mat

234:    Input Parameters:
235: +  eps - the eigenproblem solver context
236: .  A  - the matrix associated with the eigensystem
237: -  B  - the second matrix in the case of generalized eigenproblems

239:    Notes: 
240:    To specify a standard eigenproblem, use PETSC_NULL for parameter B.

242:    Level: beginner

244: .seealso: EPSSolve(), EPSGetST(), STGetOperators()
245: @*/
246: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
247: {
249:   PetscInt       m,n;


258:   /* Check for square matrices */
259:   MatGetSize(A,&m,&n);
260:   if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
261:   if (B) {
262:     MatGetSize(B,&m,&n);
263:     if (m!=n) { SETERRQ(1,"B is a non-square matrix"); }
264:   }

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

269:   /* Destroy randomly generated initial vectors */
270:   if (eps->vec_initial) {
271:     VecDestroy(eps->vec_initial);
272:     eps->vec_initial = PETSC_NULL;
273:   }
274:   if (eps->vec_initial_left) {
275:     VecDestroy(eps->vec_initial_left);
276:     eps->vec_initial_left = PETSC_NULL;
277:   }

279:   return(0);
280: }

284: /*@
285:    EPSAttachDeflationSpace - Add vectors to the basis of the deflation space.

287:    Not Collective

289:    Input Parameter:
290: +  eps   - the eigenproblem solver context
291: .  n     - number of vectors to add
292: .  ds    - set of basis vectors of the deflation space
293: -  ortho - PETSC_TRUE if basis vectors of deflation space are orthonormal

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

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

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

307:    Level: intermediate

309: .seealso: EPSRemoveDeflationSpace()
310: @*/
311: PetscErrorCode EPSAttachDeflationSpace(EPS eps,int n,Vec *ds,PetscTruth ortho)
312: {
314:   int            i;
315:   Vec            *tvec;
316: 
319:   tvec = eps->DS;
320:   if (n+eps->nds > 0) {
321:      PetscMalloc((n+eps->nds)*sizeof(Vec), &eps->DS);
322:   }
323:   if (eps->nds > 0) {
324:     for (i=0; i<eps->nds; i++) eps->DS[i] = tvec[i];
325:     PetscFree(tvec);
326:   }
327:   for (i=0; i<n; i++) {
328:     VecDuplicate(ds[i],&eps->DS[i + eps->nds]);
329:     VecCopy(ds[i],eps->DS[i + eps->nds]);
330:   }
331:   eps->nds += n;
332:   if (!ortho) eps->ds_ortho = PETSC_FALSE;
333:   eps->setupcalled = 0;
334:   return(0);
335: }

339: /*@
340:    EPSRemoveDeflationSpace - Removes the deflation space.

342:    Not Collective

344:    Input Parameter:
345: .  eps   - the eigenproblem solver context

347:    Level: intermediate

349: .seealso: EPSAttachDeflationSpace()
350: @*/
351: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
352: {
354: 
357:   if (eps->nds > 0) {
358:     VecDestroyVecs(eps->DS, eps->nds);
359:   }
360:   eps->ds_ortho = PETSC_TRUE;
361:   eps->setupcalled = 0;
362:   return(0);
363: }