Actual source code: stsles.c

slepc-3.8.0 2017-10-20
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-2017, Universitat Politecnica de Valencia, Spain

  6:    This file is part of SLEPc.
  7:    SLEPc is distributed under a 2-clause BSD license (see LICENSE).
  8:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  9: */
 10: /*
 11:    ST interface routines related to the KSP object associated with it
 12: */

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

 16: /*
 17:    This is used to set a default type for the KSP and PC objects.
 18:    It is called at STSetFromOptions (before KSPSetFromOptions)
 19:    and also at STSetUp (in case STSetFromOptions was not called).
 20: */
 21: PetscErrorCode STSetDefaultKSP(ST st)
 22: {

 28:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
 29:   if (st->ops->setdefaultksp) { (*st->ops->setdefaultksp)(st); }
 30:   return(0);
 31: }

 33: /*
 34:    This is done by all ST types except PRECOND.
 35:    The default is an LU direct solver, or GMRES+Jacobi if matmode=shell.
 36: */
 37: PetscErrorCode STSetDefaultKSP_Default(ST st)
 38: {
 40:   PC             pc;
 41:   PCType         pctype;
 42:   KSPType        ksptype;

 45:   KSPGetPC(st->ksp,&pc);
 46:   KSPGetType(st->ksp,&ksptype);
 47:   PCGetType(pc,&pctype);
 48:   if (!pctype && !ksptype) {
 49:     if (st->shift_matrix == ST_MATMODE_SHELL) {
 50:       KSPSetType(st->ksp,KSPGMRES);
 51:       PCSetType(pc,PCJACOBI);
 52:     } else {
 53:       KSPSetType(st->ksp,KSPPREONLY);
 54:       PCSetType(pc,PCLU);
 55:     }
 56:   }
 57:   return(0);
 58: }

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

 64:    Collective on ST

 66:    Input Parameters:
 67: +  st - the spectral transformation context
 68: .  k  - index of matrix to use
 69: -  x  - the vector to be multiplied

 71:    Output Parameter:
 72: .  y - the result

 74:    Level: developer

 76: .seealso: STMatMultTranspose()
 77: @*/
 78: PetscErrorCode STMatMult(ST st,PetscInt k,Vec x,Vec y)
 79: {

 87:   STCheckMatrices(st,1);
 88:   if (k<0 || k>=PetscMax(2,st->nmat)) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",st->nmat);
 89:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
 90:   VecLocked(y,3);

 92:   if (st->state!=ST_STATE_SETUP) { STSetUp(st); }
 93:   VecLockPush(x);
 94:   PetscLogEventBegin(ST_MatMult,st,x,y,0);
 95:   if (!st->T[k]) {
 96:     /* T[k]=NULL means identity matrix */
 97:     VecCopy(x,y);
 98:   } else {
 99:     MatMult(st->T[k],x,y);
100:   }
101:   PetscLogEventEnd(ST_MatMult,st,x,y,0);
102:   VecLockPop(x);
103:   return(0);
104: }

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

110:    Collective on ST

112:    Input Parameters:
113: +  st - the spectral transformation context
114: .  k  - index of matrix to use
115: -  x  - the vector to be multiplied

117:    Output Parameter:
118: .  y - the result

120:    Level: developer

122: .seealso: STMatMult()
123: @*/
124: PetscErrorCode STMatMultTranspose(ST st,PetscInt k,Vec x,Vec y)
125: {

133:   STCheckMatrices(st,1);
134:   if (k<0 || k>=PetscMax(2,st->nmat)) SETERRQ1(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %D",st->nmat);
135:   if (x == y) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
136:   VecLocked(y,3);

138:   if (st->state!=ST_STATE_SETUP) { STSetUp(st); }
139:   VecLockPush(x);
140:   PetscLogEventBegin(ST_MatMultTranspose,st,x,y,0);
141:   if (!st->T[k]) {
142:     /* T[k]=NULL means identity matrix */
143:     VecCopy(x,y);
144:   } else {
145:     MatMultTranspose(st->T[k],x,y);
146:   }
147:   PetscLogEventEnd(ST_MatMultTranspose,st,x,y,0);
148:   VecLockPop(x);
149:   return(0);
150: }

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

156:    Collective on ST

158:    Input Parameters:
159: +  st - the spectral transformation context
160: -  b  - right hand side vector

162:    Output Parameter:
163: .  x - computed solution

165:    Level: developer

167: .seealso: STMatSolveTranspose()
168: @*/
169: PetscErrorCode STMatSolve(ST st,Vec b,Vec x)
170: {
172:   PetscInt       its;
173:   PetscBool      flg;

179:   STCheckMatrices(st,1);
180:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
181:   VecLocked(x,3);

183:   if (st->state!=ST_STATE_SETUP) { STSetUp(st); }
184:   VecLockPush(b);
185:   PetscLogEventBegin(ST_MatSolve,st,b,x,0);
186:   PetscObjectTypeCompareAny((PetscObject)st,&flg,STPRECOND,STSHELL,"");
187:   if (!flg && !st->P) {
188:     /* P=NULL means identity matrix */
189:     VecCopy(b,x);
190:     return(0);
191:   }
192:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
193:   KSPSolve(st->ksp,b,x);
194:   KSPGetIterationNumber(st->ksp,&its);
195:   PetscInfo1(st,"Linear solve iterations=%D\n",its);
196:   PetscLogEventEnd(ST_MatSolve,st,b,x,0);
197:   VecLockPop(b);
198:   return(0);
199: }

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

205:    Collective on ST

207:    Input Parameters:
208: .  st - the spectral transformation context
209: .  b  - right hand side vector

211:    Output Parameter:
212: .  x - computed solution

214:    Level: developer

216: .seealso: STMatSolve()
217: @*/
218: PetscErrorCode STMatSolveTranspose(ST st,Vec b,Vec x)
219: {
221:   PetscInt       its;
222:   PetscBool      flg;

228:   STCheckMatrices(st,1);
229:   if (x == b) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
230:   VecLocked(x,3);

232:   if (st->state!=ST_STATE_SETUP) { STSetUp(st); }
233:   VecLockPush(b);
234:   PetscLogEventBegin(ST_MatSolveTranspose,st,b,x,0);
235:   PetscObjectTypeCompareAny((PetscObject)st,&flg,STPRECOND,STSHELL,"");
236:   if (!flg && !st->P) {
237:     /* P=NULL means identity matrix */
238:     VecCopy(b,x);
239:     return(0);
240:   }
241:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
242:   KSPSolveTranspose(st->ksp,b,x);
243:   KSPGetIterationNumber(st->ksp,&its);
244:   PetscInfo1(st,"Linear solve iterations=%D\n",its);
245:   PetscLogEventEnd(ST_MatSolveTranspose,st,b,x,0);
246:   VecLockPop(b);
247:   return(0);
248: }

250: /*
251:    STMatSetHermitian - Sets the Hermitian flag to the ST matrix.

253:    Input Parameters:
254: .  st - the spectral transformation context
255: .  M  - matrix
256: */
257: PetscErrorCode STMatSetHermitian(ST st,Mat M)
258: {
259: #if defined(PETSC_USE_COMPLEX)
261:   PetscBool      set,aherm,mherm;
262:   PetscInt       i;
263: #endif

266: #if defined(PETSC_USE_COMPLEX)
267:   mherm = PETSC_FALSE;
268:   for (i=0;i<st->nmat;i++) {
269:     MatIsHermitianKnown(st->A[i],&set,&aherm);
270:     if (!set) aherm = PETSC_FALSE;
271:     mherm = (mherm && aherm)? PETSC_TRUE: PETSC_FALSE;
272:     if (PetscRealPart(st->sigma)==0.0) break;
273:   }
274:   mherm = (mherm && PetscImaginaryPart(st->sigma)==0.0)? PETSC_TRUE: PETSC_FALSE;
275:   MatSetOption(M,MAT_HERMITIAN,mherm);
276: #endif
277:   return(0);
278: }

280: PetscErrorCode STCheckFactorPackage(ST st)
281: {
282:   PetscErrorCode         ierr;
283:   PC                     pc;
284:   PetscMPIInt            size;
285:   PetscBool              flg;
286:   const MatSolverPackage stype;

289:   MPI_Comm_size(PetscObjectComm((PetscObject)st),&size);
290:   if (size==1) return(0);
291:   KSPGetPC(st->ksp,&pc);
292:   PCFactorGetMatSolverPackage(pc,&stype);
293:   if (stype) {   /* currently selected PC is a factorization */
294:     PetscStrcmp(stype,MATSOLVERPETSC,&flg);
295:     if (flg) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_SUP,"You chose to solve linear systems with a factorization, but in parallel runs you need to select an external package; see the users guide for details");
296:   }
297:   return(0);
298: }

300: /*@
301:    STSetKSP - Sets the KSP object associated with the spectral
302:    transformation.

304:    Collective on ST

306:    Input Parameters:
307: +  st   - the spectral transformation context
308: -  ksp  - the linear system context

310:    Level: advanced
311: @*/
312: PetscErrorCode STSetKSP(ST st,KSP ksp)
313: {

320:   PetscObjectReference((PetscObject)ksp);
321:   KSPDestroy(&st->ksp);
322:   st->ksp = ksp;
323:   PetscLogObjectParent((PetscObject)st,(PetscObject)st->ksp);
324:   return(0);
325: }

327: /*@
328:    STGetKSP - Gets the KSP object associated with the spectral
329:    transformation.

331:    Not Collective

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

336:    Output Parameter:
337: .  ksp  - the linear system context

339:    Level: intermediate
340: @*/
341: PetscErrorCode STGetKSP(ST st,KSP* ksp)
342: {
344:   PetscBool      flg;

349:   if (!st->ksp) {
350:     KSPCreate(PetscObjectComm((PetscObject)st),&st->ksp);
351:     KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
352:     KSPAppendOptionsPrefix(st->ksp,"st_");
353:     PetscObjectIncrementTabLevel((PetscObject)st->ksp,(PetscObject)st,1);
354:     PetscLogObjectParent((PetscObject)st,(PetscObject)st->ksp);
355:     KSPSetTolerances(st->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
356:     PetscObjectTypeCompare((PetscObject)st,STPRECOND,&flg);
357:     if (!flg) {
358:       KSPSetErrorIfNotConverged(st->ksp,PETSC_TRUE);
359:     }
360:   }
361:   *ksp = st->ksp;
362:   return(0);
363: }

365: PetscErrorCode STCheckNullSpace_Default(ST st,BV V)
366: {
368:   PetscInt       nc,i,c;
369:   PetscReal      norm;
370:   Vec            *T,w,vi;
371:   Mat            A;
372:   PC             pc;
373:   MatNullSpace   nullsp;

376:   BVGetNumConstraints(V,&nc);
377:   PetscMalloc1(nc,&T);
378:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
379:   KSPGetPC(st->ksp,&pc);
380:   PCGetOperators(pc,&A,NULL);
381:   MatCreateVecs(A,NULL,&w);
382:   c = 0;
383:   for (i=0;i<nc;i++) {
384:     BVGetColumn(V,-nc+i,&vi);
385:     MatMult(A,vi,w);
386:     VecNorm(w,NORM_2,&norm);
387:     if (norm < 1e-8) {
388:       PetscInfo2(st,"Vector %D norm=%g\n",i,(double)norm);
389:       BVCreateVec(V,T+c);
390:       VecCopy(vi,T[c]);
391:       c++;
392:     }
393:     BVRestoreColumn(V,-nc+i,&vi);
394:   }
395:   VecDestroy(&w);
396:   if (c>0) {
397:     MatNullSpaceCreate(PetscObjectComm((PetscObject)st),PETSC_FALSE,c,T,&nullsp);
398:     MatSetNullSpace(A,nullsp);
399:     MatNullSpaceDestroy(&nullsp);
400:     VecDestroyVecs(c,&T);
401:   } else {
402:     PetscFree(T);
403:   }
404:   return(0);
405: }

407: /*@
408:    STCheckNullSpace - Given a basis vectors object, this function tests each
409:    of its constraint vectors to be a nullspace vector of the coefficient
410:    matrix of the associated KSP object. All these nullspace vectors are passed
411:    to the KSP object.

413:    Collective on ST

415:    Input Parameters:
416: +  st - the spectral transformation context
417: -  V  - basis vectors to be checked

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

423:    Level: developer

425: .seealso: EPSSetDeflationSpace()
426: @*/
427: PetscErrorCode STCheckNullSpace(ST st,BV V)
428: {
430:   PetscInt       nc;

437:   if (!st->state) SETERRQ(PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Must call STSolve() first");

439:   BVGetNumConstraints(V,&nc);
440:   if (nc && st->ops->checknullspace) {
441:     (*st->ops->checknullspace)(st,V);
442:   }
443:   return(0);
444: }