Actual source code: stsles.c

slepc-3.19.2 2023-09-05
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: {
 23:   PetscFunctionBegin;
 26:   if (!st->ksp) PetscCall(STGetKSP(st,&st->ksp));
 27:   PetscTryTypeMethod(st,setdefaultksp);
 28:   PetscFunctionReturn(PETSC_SUCCESS);
 29: }

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

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

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

 65:    Neighbor-wise Collective

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

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

 75:    Level: developer

 77: .seealso: STMatMultTranspose()
 78: @*/
 79: PetscErrorCode STMatMult(ST st,PetscInt k,Vec x,Vec y)
 80: {
 81:   PetscFunctionBegin;
 86:   STCheckMatrices(st,1);
 87:   PetscCheck(k>=0 && k<PetscMax(2,st->nmat),PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,st->nmat);
 88:   PetscCheck(x!=y,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
 89:   PetscCall(VecSetErrorIfLocked(y,3));

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

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

105:    Neighbor-wise Collective

107:    Input Parameters:
108: +  st - the spectral transformation context
109: .  k  - index of matrix to use
110: -  x  - the vector to be multiplied

112:    Output Parameter:
113: .  y - the result

115:    Level: developer

117: .seealso: STMatMult()
118: @*/
119: PetscErrorCode STMatMultTranspose(ST st,PetscInt k,Vec x,Vec y)
120: {
121:   PetscFunctionBegin;
126:   STCheckMatrices(st,1);
127:   PetscCheck(k>=0 && k<PetscMax(2,st->nmat),PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_OUTOFRANGE,"k must be between 0 and %" PetscInt_FMT,st->nmat);
128:   PetscCheck(x!=y,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and y must be different vectors");
129:   PetscCall(VecSetErrorIfLocked(y,3));

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

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

145:    Collective

147:    Input Parameters:
148: +  st - the spectral transformation context
149: -  b  - right hand side vector

151:    Output Parameter:
152: .  x - computed solution

154:    Level: developer

156: .seealso: STMatSolveTranspose()
157: @*/
158: PetscErrorCode STMatSolve(ST st,Vec b,Vec x)
159: {
160:   PetscFunctionBegin;
164:   STCheckMatrices(st,1);
165:   PetscCheck(x!=b,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
166:   PetscCall(VecSetErrorIfLocked(x,3));

168:   if (st->state!=ST_STATE_SETUP) PetscCall(STSetUp(st));
169:   PetscCall(VecLockReadPush(b));
170:   PetscCall(PetscLogEventBegin(ST_MatSolve,st,b,x,0));
171:   if (!st->P) PetscCall(VecCopy(b,x)); /* P=NULL means identity matrix */
172:   else PetscCall(KSPSolve(st->ksp,b,x));
173:   PetscCall(PetscLogEventEnd(ST_MatSolve,st,b,x,0));
174:   PetscCall(VecLockReadPop(b));
175:   PetscFunctionReturn(PETSC_SUCCESS);
176: }

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

182:    Collective

184:    Input Parameters:
185: +  st - the spectral transformation context
186: -  B  - right hand side vectors

188:    Output Parameter:
189: .  X - computed solutions

191:    Level: developer

193: .seealso: STMatSolve()
194: @*/
195: PetscErrorCode STMatMatSolve(ST st,Mat B,Mat X)
196: {
197:   PetscFunctionBegin;
201:   STCheckMatrices(st,1);

203:   if (st->state!=ST_STATE_SETUP) PetscCall(STSetUp(st));
204:   PetscCall(PetscLogEventBegin(ST_MatSolve,st,B,X,0));
205:   if (!st->P) PetscCall(MatCopy(B,X,SAME_NONZERO_PATTERN)); /* P=NULL means identity matrix */
206:   else PetscCall(KSPMatSolve(st->ksp,B,X));
207:   PetscCall(PetscLogEventEnd(ST_MatSolve,st,B,X,0));
208:   PetscFunctionReturn(PETSC_SUCCESS);
209: }

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

215:    Collective

217:    Input Parameters:
218: +  st - the spectral transformation context
219: -  b  - right hand side vector

221:    Output Parameter:
222: .  x - computed solution

224:    Level: developer

226: .seealso: STMatSolve()
227: @*/
228: PetscErrorCode STMatSolveTranspose(ST st,Vec b,Vec x)
229: {
230:   PetscFunctionBegin;
234:   STCheckMatrices(st,1);
235:   PetscCheck(x!=b,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_IDN,"x and b must be different vectors");
236:   PetscCall(VecSetErrorIfLocked(x,3));

238:   if (st->state!=ST_STATE_SETUP) PetscCall(STSetUp(st));
239:   PetscCall(VecLockReadPush(b));
240:   PetscCall(PetscLogEventBegin(ST_MatSolveTranspose,st,b,x,0));
241:   if (!st->P) PetscCall(VecCopy(b,x)); /* P=NULL means identity matrix */
242:   else PetscCall(KSPSolveTranspose(st->ksp,b,x));
243:   PetscCall(PetscLogEventEnd(ST_MatSolveTranspose,st,b,x,0));
244:   PetscCall(VecLockReadPop(b));
245:   PetscFunctionReturn(PETSC_SUCCESS);
246: }

248: PetscErrorCode STCheckFactorPackage(ST st)
249: {
250:   PC             pc;
251:   PetscMPIInt    size;
252:   PetscBool      flg;
253:   MatSolverType  stype;

255:   PetscFunctionBegin;
256:   PetscCallMPI(MPI_Comm_size(PetscObjectComm((PetscObject)st),&size));
257:   if (size==1) PetscFunctionReturn(PETSC_SUCCESS);
258:   PetscCall(KSPGetPC(st->ksp,&pc));
259:   PetscCall(PCFactorGetMatSolverType(pc,&stype));
260:   if (stype) {   /* currently selected PC is a factorization */
261:     PetscCall(PetscStrcmp(stype,MATSOLVERPETSC,&flg));
262:     PetscCheck(!flg,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");
263:   }
264:   PetscFunctionReturn(PETSC_SUCCESS);
265: }

267: /*@
268:    STSetKSP - Sets the KSP object associated with the spectral
269:    transformation.

271:    Collective

273:    Input Parameters:
274: +  st   - the spectral transformation context
275: -  ksp  - the linear system context

277:    Level: advanced

279: .seealso: STGetKSP()
280: @*/
281: PetscErrorCode STSetKSP(ST st,KSP ksp)
282: {
283:   PetscFunctionBegin;
286:   PetscCheckSameComm(st,1,ksp,2);
287:   STCheckNotSeized(st,1);
288:   PetscCall(PetscObjectReference((PetscObject)ksp));
289:   PetscCall(KSPDestroy(&st->ksp));
290:   st->ksp = ksp;
291:   PetscFunctionReturn(PETSC_SUCCESS);
292: }

294: /*@
295:    STGetKSP - Gets the KSP object associated with the spectral
296:    transformation.

298:    Collective

300:    Input Parameter:
301: .  st - the spectral transformation context

303:    Output Parameter:
304: .  ksp  - the linear system context

306:    Level: intermediate

308: .seealso: STSetKSP()
309: @*/
310: PetscErrorCode STGetKSP(ST st,KSP* ksp)
311: {
312:   PetscFunctionBegin;
315:   if (!st->ksp) {
316:     PetscCall(KSPCreate(PetscObjectComm((PetscObject)st),&st->ksp));
317:     PetscCall(PetscObjectIncrementTabLevel((PetscObject)st->ksp,(PetscObject)st,1));
318:     PetscCall(KSPSetOptionsPrefix(st->ksp,((PetscObject)st)->prefix));
319:     PetscCall(KSPAppendOptionsPrefix(st->ksp,"st_"));
320:     PetscCall(PetscObjectSetOptions((PetscObject)st->ksp,((PetscObject)st)->options));
321:     PetscCall(KSPSetTolerances(st->ksp,SLEPC_DEFAULT_TOL,PETSC_DEFAULT,PETSC_DEFAULT,PETSC_DEFAULT));
322:   }
323:   *ksp = st->ksp;
324:   PetscFunctionReturn(PETSC_SUCCESS);
325: }

327: PetscErrorCode STCheckNullSpace_Default(ST st,BV V)
328: {
329:   PetscInt       nc,i,c;
330:   PetscReal      norm;
331:   Vec            *T,w,vi;
332:   Mat            A;
333:   PC             pc;
334:   MatNullSpace   nullsp;

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

367: /*@
368:    STCheckNullSpace - Tests if constraint vectors are nullspace vectors.

370:    Collective

372:    Input Parameters:
373: +  st - the spectral transformation context
374: -  V  - basis vectors to be checked

376:    Notes:
377:    Given a basis vectors object, this function tests each of its constraint
378:    vectors to be a nullspace vector of the coefficient matrix of the
379:    associated KSP object. All these nullspace vectors are passed to the KSP
380:    object.

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

385:    Level: developer

387: .seealso: EPSSetDeflationSpace()
388: @*/
389: PetscErrorCode STCheckNullSpace(ST st,BV V)
390: {
391:   PetscInt       nc;

393:   PetscFunctionBegin;
397:   PetscCheckSameComm(st,1,V,2);
398:   PetscCheck(st->state,PetscObjectComm((PetscObject)st),PETSC_ERR_ARG_WRONGSTATE,"Must call STSetUp() first");

400:   PetscCall(BVGetNumConstraints(V,&nc));
401:   if (nc) PetscTryTypeMethod(st,checknullspace,V);
402:   PetscFunctionReturn(PETSC_SUCCESS);
403: }