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:   PetscTruth     isCayley;
 57: 

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

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

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

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

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

114:   PetscTypeCompare((PetscObject)eps->OP,STCAYLEY,&isCayley);
115:   if (isCayley && eps->problem_type == EPS_PGNHEP) {
116:     SETERRQ(PETSC_ERR_SUP,"Cayley spectral transformation is not compatible with PGNHEP");
117:   }

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

133:   STCheckNullSpace(eps->OP,eps->nds,eps->DS);
134: 
135:   PetscLogEventEnd(EPS_SetUp,eps,0,0,0);
136:   eps->setupcalled = 1;
137:   return(0);
138: }

142: /*@
143:    EPSSetInitialVector - Sets the initial vector from which the 
144:    eigensolver starts to iterate.

146:    Collective on EPS and Vec

148:    Input Parameters:
149: +  eps - the eigensolver context
150: -  vec - the vector

152:    Level: intermediate

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

156: @*/
157: PetscErrorCode EPSSetInitialVector(EPS eps,Vec vec)
158: {
160: 
165:   PetscObjectReference((PetscObject)vec);
166:   if (eps->vec_initial) {
167:     VecDestroy(eps->vec_initial);
168:   }
169:   eps->vec_initial = vec;
170:   return(0);
171: }

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

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

182:    Input Parameter:
183: .  eps - the eigensolver context

185:    Output Parameter:
186: .  vec - the vector

188:    Level: intermediate

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

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

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

208:    Collective on EPS and Vec

210:    Input Parameters:
211: +  eps - the eigensolver context
212: -  vec - the vector

214:    Level: intermediate

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

218: @*/
219: PetscErrorCode EPSSetLeftInitialVector(EPS eps,Vec vec)
220: {
222: 
227:   PetscObjectReference((PetscObject)vec);
228:   if (eps->vec_initial_left) {
229:     VecDestroy(eps->vec_initial_left);
230:   }
231:   eps->vec_initial_left = vec;
232:   return(0);
233: }

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

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

244:    Input Parameter:
245: .  eps - the eigensolver context

247:    Output Parameter:
248: .  vec - the vector

250:    Level: intermediate

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

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

266: /*@
267:    EPSSetOperators - Sets the matrices associated with the eigenvalue problem.

269:    Collective on EPS and Mat

271:    Input Parameters:
272: +  eps - the eigenproblem solver context
273: .  A  - the matrix associated with the eigensystem
274: -  B  - the second matrix in the case of generalized eigenproblems

276:    Notes: 
277:    To specify a standard eigenproblem, use PETSC_NULL for parameter B.

279:    Level: beginner

281: .seealso: EPSSolve(), EPSGetST(), STGetOperators()
282: @*/
283: PetscErrorCode EPSSetOperators(EPS eps,Mat A,Mat B)
284: {
286:   PetscInt       m,n;


295:   /* Check for square matrices */
296:   MatGetSize(A,&m,&n);
297:   if (m!=n) { SETERRQ(1,"A is a non-square matrix"); }
298:   if (B) {
299:     MatGetSize(B,&m,&n);
300:     if (m!=n) { SETERRQ(1,"B is a non-square matrix"); }
301:   }

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

306:   /* Destroy randomly generated initial vectors */
307:   if (eps->vec_initial) {
308:     VecDestroy(eps->vec_initial);
309:     eps->vec_initial = PETSC_NULL;
310:   }
311:   if (eps->vec_initial_left) {
312:     VecDestroy(eps->vec_initial_left);
313:     eps->vec_initial_left = PETSC_NULL;
314:   }

316:   return(0);
317: }

321: /*@
322:    EPSGetOperators - Gets the matrices associated with the eigensystem.

324:    Collective on EPS and Mat

326:    Input Parameter:
327: .  eps - the EPS context

329:    Output Parameters:
330: +  A  - the matrix associated with the eigensystem
331: -  B  - the second matrix in the case of generalized eigenproblems

333:    Level: intermediate

335: .seealso: EPSSolve(), EPSGetST(), STGetOperators(), STSetOperators()
336: @*/
337: PetscErrorCode EPSGetOperators(EPS eps, Mat *A, Mat *B)
338: {
340:   ST             st;

346:   EPSGetST(eps,&st);
347:   STGetOperators(st,A,B);
348:   return(0);
349: }

353: /*@
354:    EPSAttachDeflationSpace - Add vectors to the basis of the deflation space.

356:    Not Collective

358:    Input Parameter:
359: +  eps   - the eigenproblem solver context
360: .  n     - number of vectors to add
361: .  ds    - set of basis vectors of the deflation space
362: -  ortho - PETSC_TRUE if basis vectors of deflation space are orthonormal

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

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

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

376:    Level: intermediate

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

408: /*@
409:    EPSRemoveDeflationSpace - Removes the deflation space.

411:    Not Collective

413:    Input Parameter:
414: .  eps   - the eigenproblem solver context

416:    Level: intermediate

418: .seealso: EPSAttachDeflationSpace()
419: @*/
420: PetscErrorCode EPSRemoveDeflationSpace(EPS eps)
421: {
423: 
426:   if (eps->nds > 0) {
427:     VecDestroyVecs(eps->DS, eps->nds);
428:   }
429:   eps->ds_ortho = PETSC_TRUE;
430:   eps->setupcalled = 0;
431:   return(0);
432: }