Actual source code: stsles.c

slepc-3.5.2 2014-10-10
Report Typos and Errors
  1: /*
  2:     The ST (spectral transformation) interface routines related to the
  3:     KSP object associated to it.

  5:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  6:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  7:    Copyright (c) 2002-2014, Universitat Politecnica de Valencia, Spain

  9:    This file is part of SLEPc.

 11:    SLEPc is free software: you can redistribute it and/or modify it under  the
 12:    terms of version 3 of the GNU Lesser General Public License as published by
 13:    the Free Software Foundation.

 15:    SLEPc  is  distributed in the hope that it will be useful, but WITHOUT  ANY
 16:    WARRANTY;  without even the implied warranty of MERCHANTABILITY or  FITNESS
 17:    FOR  A  PARTICULAR PURPOSE. See the GNU Lesser General Public  License  for
 18:    more details.

 20:    You  should have received a copy of the GNU Lesser General  Public  License
 21:    along with SLEPc. If not, see <http://www.gnu.org/licenses/>.
 22:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 23: */

 25: #include <slepc-private/stimpl.h>            /*I "slepcst.h" I*/

 29: /*@
 30:    STMatMult - Computes the matrix-vector product y = T[k] x, where T[k] is
 31:    the k-th matrix of the spectral transformation.

 33:    Collective on ST

 35:    Input Parameters:
 36: +  st - the spectral transformation context
 37: .  k  - index of matrix to use
 38: -  x  - the vector to be multiplied

 40:    Output Parameter:
 41: .  y - the result

 43:    Level: developer

 45: .seealso: STMatMultTranspose()
 46: @*/
 47: PetscErrorCode STMatMult(ST st,PetscInt k,Vec x,Vec y)
 48: {

 56:   STCheckMatrices(st,1);
 57:   if (k<0 || k>=PetscMax(2,st->nmat)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",st->nmat);
 58:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");

 60:   if (!st->setupcalled) { STSetUp(st); }
 61:   PetscLogEventBegin(ST_MatMult,st,x,y,0);
 62:   if (!st->T[k]) {
 63:     /* T[k]=NULL means identity matrix */
 64:     VecCopy(x,y);
 65:   } else {
 66:     MatMult(st->T[k],x,y);
 67:   }
 68:   PetscLogEventEnd(ST_MatMult,st,x,y,0);
 69:   return(0);
 70: }

 74: /*@
 75:    STMatMultTranspose - Computes the matrix-vector product y = T[k]' x, where T[k] is
 76:    the k-th matrix of the spectral transformation.

 78:    Collective on ST

 80:    Input Parameters:
 81: +  st - the spectral transformation context
 82: .  k  - index of matrix to use
 83: -  x  - the vector to be multiplied

 85:    Output Parameter:
 86: .  y - the result

 88:    Level: developer

 90: .seealso: STMatMult()
 91: @*/
 92: PetscErrorCode STMatMultTranspose(ST st,PetscInt k,Vec x,Vec y)
 93: {

101:   STCheckMatrices(st,1);
102:   if (k<0 || k>=PetscMax(2,st->nmat)) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %d",st->nmat);
103:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");

105:   if (!st->setupcalled) { STSetUp(st); }
106:   PetscLogEventBegin(ST_MatMultTranspose,st,x,y,0);
107:   if (!st->T[k]) {
108:     /* T[k]=NULL means identity matrix */
109:     VecCopy(x,y);
110:   } else {
111:     MatMultTranspose(st->T[k],x,y);
112:   }
113:   PetscLogEventEnd(ST_MatMultTranspose,st,x,y,0);
114:   return(0);
115: }

119: /*@
120:    STMatSolve - Solves P x = b, where P is the preconditioner matrix of
121:    the spectral transformation, using a KSP object stored internally.

123:    Collective on ST

125:    Input Parameters:
126: +  st - the spectral transformation context
127: -  b  - right hand side vector

129:    Output Parameter:
130: .  x - computed solution

132:    Level: developer

134: .seealso: STMatSolveTranspose()
135: @*/
136: PetscErrorCode STMatSolve(ST st,Vec b,Vec x)
137: {
138:   PetscErrorCode     ierr;
139:   PetscInt           its;
140:   PetscBool          flg;
141:   KSPConvergedReason reason;

147:   STCheckMatrices(st,1);
148:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and b must be different vectors");

150:   if (!st->setupcalled) { STSetUp(st); }
151:   PetscLogEventBegin(ST_MatSolve,st,b,x,0);
152:   PetscObjectTypeCompareAny((PetscObject)st,&flg,STPRECOND,STSHELL,"");
153:   if (!flg && !st->P) {
154:     /* P=NULL means identity matrix */
155:     VecCopy(b,x);
156:     return(0);
157:   }
158:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
159:   KSPSolve(st->ksp,b,x);
160:   KSPGetConvergedReason(st->ksp,&reason);
161:   if (reason<0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
162:   KSPGetIterationNumber(st->ksp,&its);
163:   st->linearits += its;
164:   PetscInfo1(st,"Linear solve iterations=%D\n",its);
165:   PetscLogEventEnd(ST_MatSolve,st,b,x,0);
166:   return(0);
167: }

171: /*@
172:    STMatSolveTranspose - Solves P' x = b, where P is the preconditioner matrix of
173:    the spectral transformation, using a KSP object stored internally.

175:    Collective on ST

177:    Input Parameters:
178: .  st - the spectral transformation context
179: .  b  - right hand side vector

181:    Output Parameter:
182: .  x - computed solution

184:    Level: developer

186: .seealso: STMatSolve()
187: @*/
188: PetscErrorCode STMatSolveTranspose(ST st,Vec b,Vec x)
189: {
190:   PetscErrorCode     ierr;
191:   PetscInt           its;
192:   PetscBool          flg;
193:   KSPConvergedReason reason;

199:   STCheckMatrices(st,1);
200:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and b must be different vectors");

202:   if (!st->setupcalled) { STSetUp(st); }
203:   PetscLogEventBegin(ST_MatSolveTranspose,st,b,x,0);
204:   PetscObjectTypeCompareAny((PetscObject)st,&flg,STPRECOND,STSHELL,"");
205:   if (!flg && !st->P) {
206:     /* P=NULL means identity matrix */
207:     VecCopy(b,x);
208:     return(0);
209:   }
210:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
211:   KSPSolveTranspose(st->ksp,b,x);
212:   KSPGetConvergedReason(st->ksp,&reason);
213:   if (reason<0) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_NOT_CONVERGED,"KSP did not converge (reason=%s)",KSPConvergedReasons[reason]);
214:   KSPGetIterationNumber(st->ksp,&its);
215:   st->linearits += its;
216:   PetscInfo1(st,"Linear solve iterations=%D\n",its);
217:   PetscLogEventEnd(ST_MatSolveTranspose,st,b,x,0);
218:   return(0);
219: }

223: /*
224:    STMatSetHermitian - Sets the Hermitian flag to the ST matrix.

226:    Input Parameters:
227: .  st - the spectral transformation context
228: .  M  - matrix
229: */
230: PetscErrorCode STMatSetHermitian(ST st,Mat M)
231: {
232: #if defined(PETSC_USE_COMPLEX)
234:   PetscBool      set,aherm,mherm;
235:   PetscInt       i;
236: #endif

239: #if defined(PETSC_USE_COMPLEX)
240:   mherm = PETSC_FALSE;
241:   for (i=0;i<st->nmat;i++) {
242:     MatIsHermitianKnown(st->A[i],&set,&aherm);
243:     if (!set) aherm = PETSC_FALSE;
244:     mherm = (mherm && aherm)? PETSC_TRUE: PETSC_FALSE;
245:     if (PetscRealPart(st->sigma)==0.0) break;
246:   }
247:   mherm = (mherm && PetscImaginaryPart(st->sigma)==0.0)? PETSC_TRUE: PETSC_FALSE;
248:   MatSetOption(M,MAT_HERMITIAN,mherm);
249: #endif
250:   return(0);
251: }

255: /*@
256:    STSetKSP - Sets the KSP object associated with the spectral
257:    transformation.

259:    Collective on ST

261:    Input Parameters:
262: +  st   - the spectral transformation context
263: -  ksp  - the linear system context

265:    Level: advanced

267: @*/
268: PetscErrorCode STSetKSP(ST st,KSP ksp)
269: {

276:   PetscObjectReference((PetscObject)ksp);
277:   KSPDestroy(&st->ksp);
278:   st->ksp = ksp;
279:   PetscLogObjectParent((PetscObject)st,(PetscObject)st->ksp);
280:   return(0);
281: }

285: /*@
286:    STGetKSP - Gets the KSP object associated with the spectral
287:    transformation.

289:    Not Collective

291:    Input Parameter:
292: .  st - the spectral transformation context

294:    Output Parameter:
295: .  ksp  - the linear system context

297:    Notes:
298:    On output, the value of ksp can be NULL if the combination of
299:    eigenproblem type and selected transformation does not require to
300:    solve a linear system of equations.

302:    Level: intermediate

304: @*/
305: PetscErrorCode STGetKSP(ST st,KSP* ksp)
306: {

312:   if (!st->ksp) {
313:     KSPCreate(PetscObjectComm((PetscObject)st),&st->ksp);
314:     KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
315:     KSPAppendOptionsPrefix(st->ksp,"st_");
316:     PetscObjectIncrementTabLevel((PetscObject)st->ksp,(PetscObject)st,1);
317:     PetscLogObjectParent((PetscObject)st,(PetscObject)st->ksp);
318:     KSPSetTolerances(st->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
319:   }
320:   *ksp = st->ksp;
321:   return(0);
322: }

326: /*@
327:    STGetOperationCounters - Gets the total number of operator applications
328:    and linear solver iterations used by the ST object.

330:    Not Collective

332:    Input Parameter:
333: .  st - the spectral transformation context

335:    Output Parameter:
336: +  ops  - number of operator applications
337: -  lits - number of linear solver iterations

339:    Notes:
340:    Any output parameter may be NULL on input if not needed.

342:    Level: intermediate

344: .seealso: STResetOperationCounters()
345: @*/
346: PetscErrorCode STGetOperationCounters(ST st,PetscInt* ops,PetscInt* lits)
347: {
350:   if (ops)  *ops  = st->applys;
351:   if (lits) *lits = st->linearits;
352:   return(0);
353: }

357: /*@
358:    STResetOperationCounters - Resets the counters for operator applications,
359:    inner product operations and total number of linear iterations used by
360:    the ST object.

362:    Logically Collective on ST

364:    Input Parameter:
365: .  st - the spectral transformation context

367:    Level: intermediate

369: .seealso: STGetOperationCounters()
370: @*/
371: PetscErrorCode STResetOperationCounters(ST st)
372: {
375:   st->linearits = 0;
376:   st->applys    = 0;
377:   return(0);
378: }

382: PetscErrorCode STCheckNullSpace_Default(ST st,BV V)
383: {
385:   PetscInt       nc,i,c;
386:   PetscReal      norm;
387:   Vec            *T,w,vi;
388:   Mat            A;
389:   PC             pc;
390:   MatNullSpace   nullsp;

393:   BVGetNumConstraints(V,&nc);
394:   PetscMalloc1(nc,&T);
395:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
396:   KSPGetPC(st->ksp,&pc);
397:   PCGetOperators(pc,&A,NULL);
398:   MatGetVecs(A,NULL,&w);
399:   c = 0;
400:   for (i=0;i<nc;i++) {
401:     BVGetColumn(V,-nc+i,&vi);
402:     MatMult(A,vi,w);
403:     VecNorm(w,NORM_2,&norm);
404:     if (norm < 1e-8) {
405:       PetscInfo2(st,"Vector %D norm=%g\n",i,(double)norm);
406:       BVGetVec(V,T+c);
407:       VecCopy(vi,T[c]);
408:       c++;
409:     }
410:     BVRestoreColumn(V,-nc+i,&vi);
411:   }
412:   VecDestroy(&w);
413:   if (c>0) {
414:     MatNullSpaceCreate(PetscObjectComm((PetscObject)st),PETSC_FALSE,c,T,&nullsp);
415:     KSPSetNullSpace(st->ksp,nullsp);
416:     MatNullSpaceDestroy(&nullsp);
417:     VecDestroyVecs(c,&T);
418:   } else {
419:     PetscFree(T);
420:   }
421:   return(0);
422: }

426: /*@
427:    STCheckNullSpace - Given a basis vectors object, this function tests each
428:    of its constraint vectors to be a nullspace vector of the coefficient
429:    matrix of the associated KSP object. All these nullspace vectors are passed
430:    to the KSP object.

432:    Collective on ST

434:    Input Parameters:
435: +  st - the spectral transformation context
436: -  V  - basis vectors to be checked

438:    Note:
439:    This function allows to handle singular pencils and to solve some problems
440:    in which the nullspace is important (see the users guide for details).

442:    Level: developer

444: .seealso: EPSSetDeflationSpace()
445: @*/
446: PetscErrorCode STCheckNullSpace(ST st,BV V)
447: {
449:   PetscInt       nc;


457:   BVGetNumConstraints(V,&nc);
458:   if (nc && st->ops->checknullspace) {
459:     (*st->ops->checknullspace)(st,V);
460:   }
461:   return(0);
462: }