Actual source code: stsles.c

slepc-3.17.1 2022-04-11
Report Typos and Errors
  1: /*
  2:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
  3:    SLEPc - Scalable Library for Eigenvalue Problem Computations
  4:    Copyright (c) 2002-, 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>

 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: {
 25:   if (!st->ksp) STGetKSP(st,&st->ksp);
 26:   if (st->ops->setdefaultksp) (*st->ops->setdefaultksp)(st);
 27:   PetscFunctionReturn(0);
 28: }

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

 40:   KSPGetPC(st->ksp,&pc);
 41:   KSPGetType(st->ksp,&ksptype);
 42:   PCGetType(pc,&pctype);
 43:   if (!pctype && !ksptype) {
 44:     if (st->Pmat || st->Psplit) {
 45:       KSPSetType(st->ksp,KSPBCGS);
 46:       PCSetType(pc,PCBJACOBI);
 47:     } else if (st->matmode == ST_MATMODE_SHELL) {
 48:       KSPSetType(st->ksp,KSPGMRES);
 49:       PCSetType(pc,PCJACOBI);
 50:     } else {
 51:       KSPSetType(st->ksp,KSPPREONLY);
 52:       PCSetType(pc,st->asymm?PCCHOLESKY:PCLU);
 53:     }
 54:   }
 55:   KSPSetErrorIfNotConverged(st->ksp,PETSC_TRUE);
 56:   PetscFunctionReturn(0);
 57: }

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

 63:    Neighbor-wise Collective on st

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

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

 73:    Level: developer

 75: .seealso: STMatMultTranspose()
 76: @*/
 77: PetscErrorCode STMatMult(ST st,PetscInt k,Vec x,Vec y)
 78: {
 83:   STCheckMatrices(st,1);
 86:   VecSetErrorIfLocked(y,3);

 88:   if (st->state!=ST_STATE_SETUP) STSetUp(st);
 89:   VecLockReadPush(x);
 90:   PetscLogEventBegin(ST_MatMult,st,x,y,0);
 91:   if (!st->T[k]) VecCopy(x,y); /* T[k]=NULL means identity matrix */
 92:   else MatMult(st->T[k],x,y);
 93:   PetscLogEventEnd(ST_MatMult,st,x,y,0);
 94:   VecLockReadPop(x);
 95:   PetscFunctionReturn(0);
 96: }

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

102:    Neighbor-wise Collective on st

104:    Input Parameters:
105: +  st - the spectral transformation context
106: .  k  - index of matrix to use
107: -  x  - the vector to be multiplied

109:    Output Parameter:
110: .  y - the result

112:    Level: developer

114: .seealso: STMatMult()
115: @*/
116: PetscErrorCode STMatMultTranspose(ST st,PetscInt k,Vec x,Vec y)
117: {
122:   STCheckMatrices(st,1);
125:   VecSetErrorIfLocked(y,3);

127:   if (st->state!=ST_STATE_SETUP) STSetUp(st);
128:   VecLockReadPush(x);
129:   PetscLogEventBegin(ST_MatMultTranspose,st,x,y,0);
130:   if (!st->T[k]) VecCopy(x,y); /* T[k]=NULL means identity matrix */
131:   else MatMultTranspose(st->T[k],x,y);
132:   PetscLogEventEnd(ST_MatMultTranspose,st,x,y,0);
133:   VecLockReadPop(x);
134:   PetscFunctionReturn(0);
135: }

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

141:    Collective on st

143:    Input Parameters:
144: +  st - the spectral transformation context
145: -  b  - right hand side vector

147:    Output Parameter:
148: .  x - computed solution

150:    Level: developer

152: .seealso: STMatSolveTranspose()
153: @*/
154: PetscErrorCode STMatSolve(ST st,Vec b,Vec x)
155: {
159:   STCheckMatrices(st,1);
161:   VecSetErrorIfLocked(x,3);

163:   if (st->state!=ST_STATE_SETUP) STSetUp(st);
164:   VecLockReadPush(b);
165:   PetscLogEventBegin(ST_MatSolve,st,b,x,0);
166:   if (!st->P) VecCopy(b,x); /* P=NULL means identity matrix */
167:   else KSPSolve(st->ksp,b,x);
168:   PetscLogEventEnd(ST_MatSolve,st,b,x,0);
169:   VecLockReadPop(b);
170:   PetscFunctionReturn(0);
171: }

173: /*@
174:    STMatMatSolve - Solves P X = B, where P is the preconditioner matrix of
175:    the spectral transformation, using a KSP object stored internally.

177:    Collective on st

179:    Input Parameters:
180: +  st - the spectral transformation context
181: -  B  - right hand side vectors

183:    Output Parameter:
184: .  X - computed solutions

186:    Level: developer

188: .seealso: STMatSolve()
189: @*/
190: PetscErrorCode STMatMatSolve(ST st,Mat B,Mat X)
191: {
195:   STCheckMatrices(st,1);

197:   if (st->state!=ST_STATE_SETUP) STSetUp(st);
198:   PetscLogEventBegin(ST_MatSolve,st,B,X,0);
199:   if (!st->P) MatCopy(B,X,SAME_NONZERO_PATTERN); /* P=NULL means identity matrix */
200:   else KSPMatSolve(st->ksp,B,X);
201:   PetscLogEventEnd(ST_MatSolve,st,B,X,0);
202:   PetscFunctionReturn(0);
203: }

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

209:    Collective on st

211:    Input Parameters:
212: +  st - the spectral transformation context
213: -  b  - right hand side vector

215:    Output Parameter:
216: .  x - computed solution

218:    Level: developer

220: .seealso: STMatSolve()
221: @*/
222: PetscErrorCode STMatSolveTranspose(ST st,Vec b,Vec x)
223: {
227:   STCheckMatrices(st,1);
229:   VecSetErrorIfLocked(x,3);

231:   if (st->state!=ST_STATE_SETUP) STSetUp(st);
232:   VecLockReadPush(b);
233:   PetscLogEventBegin(ST_MatSolveTranspose,st,b,x,0);
234:   if (!st->P) VecCopy(b,x); /* P=NULL means identity matrix */
235:   else KSPSolveTranspose(st->ksp,b,x);
236:   PetscLogEventEnd(ST_MatSolveTranspose,st,b,x,0);
237:   VecLockReadPop(b);
238:   PetscFunctionReturn(0);
239: }

241: PetscErrorCode STCheckFactorPackage(ST st)
242: {
243:   PC             pc;
244:   PetscMPIInt    size;
245:   PetscBool      flg;
246:   MatSolverType  stype;

248:   MPI_Comm_size(PetscObjectComm((PetscObject)st),&size);
249:   if (size==1) PetscFunctionReturn(0);
250:   KSPGetPC(st->ksp,&pc);
251:   PCFactorGetMatSolverType(pc,&stype);
252:   if (stype) {   /* currently selected PC is a factorization */
253:     PetscStrcmp(stype,MATSOLVERPETSC,&flg);
255:   }
256:   PetscFunctionReturn(0);
257: }

259: /*@
260:    STSetKSP - Sets the KSP object associated with the spectral
261:    transformation.

263:    Collective on st

265:    Input Parameters:
266: +  st   - the spectral transformation context
267: -  ksp  - the linear system context

269:    Level: advanced

271: .seealso: STGetKSP()
272: @*/
273: PetscErrorCode STSetKSP(ST st,KSP ksp)
274: {
278:   STCheckNotSeized(st,1);
279:   PetscObjectReference((PetscObject)ksp);
280:   KSPDestroy(&st->ksp);
281:   st->ksp = ksp;
282:   PetscLogObjectParent((PetscObject)st,(PetscObject)st->ksp);
283:   PetscFunctionReturn(0);
284: }

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

290:    Not Collective

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

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

298:    Level: intermediate

300: .seealso: STSetKSP()
301: @*/
302: PetscErrorCode STGetKSP(ST st,KSP* ksp)
303: {
306:   if (!st->ksp) {
307:     KSPCreate(PetscObjectComm((PetscObject)st),&st->ksp);
308:     PetscObjectIncrementTabLevel((PetscObject)st->ksp,(PetscObject)st,1);
309:     KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix);
310:     KSPAppendOptionsPrefix(st->ksp,"st_");
311:     PetscLogObjectParent((PetscObject)st,(PetscObject)st->ksp);
312:     PetscObjectSetOptions((PetscObject)st->ksp,((PetscObject)st)->options);
313:     KSPSetTolerances(st->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT);
314:   }
315:   *ksp = st->ksp;
316:   PetscFunctionReturn(0);
317: }

319: PetscErrorCode STCheckNullSpace_Default(ST st,BV V)
320: {
321:   PetscInt       nc,i,c;
322:   PetscReal      norm;
323:   Vec            *T,w,vi;
324:   Mat            A;
325:   PC             pc;
326:   MatNullSpace   nullsp;

328:   BVGetNumConstraints(V,&nc);
329:   PetscMalloc1(nc,&T);
330:   if (!st->ksp) STGetKSP(st,&st->ksp);
331:   KSPGetPC(st->ksp,&pc);
332:   PCGetOperators(pc,&A,NULL);
333:   MatCreateVecs(A,NULL,&w);
334:   c = 0;
335:   for (i=0;i<nc;i++) {
336:     BVGetColumn(V,-nc+i,&vi);
337:     MatMult(A,vi,w);
338:     VecNorm(w,NORM_2,&norm);
339:     if (norm < 10.0*PETSC_SQRT_MACHINE_EPSILON) {
340:       PetscInfo(st,"Vector %" PetscInt_FMT " included in the nullspace of OP, norm=%g\n",i,(double)norm);
341:       BVCreateVec(V,T+c);
342:       VecCopy(vi,T[c]);
343:       VecNormalize(T[c],NULL);
344:       c++;
345:     } else PetscInfo(st,"Vector %" PetscInt_FMT " discarded as possible nullspace of OP, norm=%g\n",i,(double)norm);
346:     BVRestoreColumn(V,-nc+i,&vi);
347:   }
348:   VecDestroy(&w);
349:   if (c>0) {
350:     MatNullSpaceCreate(PetscObjectComm((PetscObject)st),PETSC_FALSE,c,T,&nullsp);
351:     MatSetNullSpace(A,nullsp);
352:     MatNullSpaceDestroy(&nullsp);
353:     VecDestroyVecs(c,&T);
354:   } else PetscFree(T);
355:   PetscFunctionReturn(0);
356: }

358: /*@
359:    STCheckNullSpace - Tests if constraint vectors are nullspace vectors.

361:    Collective on st

363:    Input Parameters:
364: +  st - the spectral transformation context
365: -  V  - basis vectors to be checked

367:    Notes:
368:    Given a basis vectors object, this function tests each of its constraint
369:    vectors to be a nullspace vector of the coefficient matrix of the
370:    associated KSP object. All these nullspace vectors are passed to the KSP
371:    object.

373:    This function allows handling singular pencils and solving some problems
374:    in which the nullspace is important (see the users guide for details).

376:    Level: developer

378: .seealso: EPSSetDeflationSpace()
379: @*/
380: PetscErrorCode STCheckNullSpace(ST st,BV V)
381: {
382:   PetscInt       nc;


390:   BVGetNumConstraints(V,&nc);
391:   if (nc && st->ops->checknullspace) (*st->ops->checknullspace)(st,V);
392:   PetscFunctionReturn(0);
393: }