Actual source code: precond.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:    Implements the ST class for preconditioned eigenvalue methods
 12: */

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

 16: typedef struct {
 17:   PetscBool setmat;
 18: } ST_PRECOND;

 20: static PetscErrorCode STSetDefaultKSP_Precond(ST st)
 21: {
 23:   PC             pc;
 24:   PCType         pctype;
 25:   Mat            P;
 26:   PetscBool      t0,t1;

 29:   KSPGetPC(st->ksp,&pc);
 30:   PetscObjectGetType((PetscObject)pc,&pctype);
 31:   STPrecondGetMatForPC(st,&P);
 32:   if (!pctype && st->A && st->A[0]) {
 33:     if (P || st->shift_matrix == ST_MATMODE_SHELL) {
 34:       PCSetType(pc,PCJACOBI);
 35:     } else {
 36:       MatHasOperation(st->A[0],MATOP_DUPLICATE,&t0);
 37:       if (st->nmat>1) {
 38:         MatHasOperation(st->A[0],MATOP_AXPY,&t1);
 39:       } else t1 = PETSC_TRUE;
 40:       PCSetType(pc,(t0 && t1)?PCJACOBI:PCNONE);
 41:     }
 42:   }
 43:   return(0);
 44: }

 46: PetscErrorCode STSetUp_Precond(ST st)
 47: {
 48:   Mat            P;
 49:   PC             pc;
 50:   PetscBool      t0,setmat,destroyP=PETSC_FALSE,builtP;

 54:   /* if the user did not set the shift, use the target value */
 55:   if (!st->sigma_set) st->sigma = st->defsigma;

 57:   /* If either pc is none and no matrix has to be set, or pc is shell , exit */
 58:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
 59:   KSPGetPC(st->ksp,&pc);
 60:   PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&t0);
 61:   if (t0) return(0);
 62:   PetscObjectTypeCompare((PetscObject)pc,PCNONE,&t0);
 63:   STPrecondGetKSPHasMat(st,&setmat);
 64:   if (t0 && !setmat) return(0);

 66:   /* Check if a user matrix is set */
 67:   STPrecondGetMatForPC(st,&P);

 69:   /* If not, create A - shift*B */
 70:   if (P) {
 71:     builtP = PETSC_FALSE;
 72:     destroyP = PETSC_TRUE;
 73:     PetscObjectReference((PetscObject)P);
 74:   } else {
 75:     builtP = PETSC_TRUE;

 77:     if (!(PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) && st->nmat>1) {
 78:       P = st->A[1];
 79:       destroyP = PETSC_FALSE;
 80:     } else if (st->sigma == 0.0) {
 81:       P = st->A[0];
 82:       destroyP = PETSC_FALSE;
 83:     } else if (PetscAbsScalar(st->sigma) < PETSC_MAX_REAL && st->shift_matrix != ST_MATMODE_SHELL) {
 84:       if (st->shift_matrix == ST_MATMODE_INPLACE) {
 85:         P = st->A[0];
 86:         destroyP = PETSC_FALSE;
 87:       } else {
 88:         MatDuplicate(st->A[0],MAT_COPY_VALUES,&P);
 89:         destroyP = PETSC_TRUE;
 90:       }
 91:       if (st->nmat>1) {
 92:         MatAXPY(P,-st->sigma,st->A[1],st->str);
 93:       } else {
 94:         MatShift(P,-st->sigma);
 95:       }
 96:       /* TODO: in case of ST_MATMODE_INPLACE should keep the Hermitian flag of st->A and restore at the end */
 97:       STMatSetHermitian(st,P);
 98:     } else builtP = PETSC_FALSE;
 99:   }

101:   /* If P was not possible to obtain, set pc to PCNONE */
102:   if (!P) {
103:     PCSetType(pc,PCNONE);

105:     /* If some matrix has to be set to ksp, a shell matrix is created */
106:     if (setmat) {
107:       STMatShellCreate(st,-st->sigma,0,NULL,NULL,&P);
108:       STMatSetHermitian(st,P);
109:       destroyP = PETSC_TRUE;
110:     }
111:   }

113:   KSPSetOperators(st->ksp,setmat?P:NULL,P);

115:   if (destroyP) {
116:     MatDestroy(&P);
117:   } else if (st->shift_matrix == ST_MATMODE_INPLACE && builtP) {
118:     if (st->sigma != 0.0 && PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) {
119:       if (st->nmat>1) {
120:         MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
121:       } else {
122:         MatShift(st->A[0],st->sigma);
123:       }
124:     }
125:   }
126:   return(0);
127: }

129: PetscErrorCode STSetShift_Precond(ST st,PetscScalar newshift)
130: {

134:   /* Nothing to be done if STSetUp has not been called yet */
135:   if (!st->state) return(0);
136:   st->sigma = newshift;
137:   if (st->shift_matrix != ST_MATMODE_SHELL) {
138:     STSetUp_Precond(st);
139:   }
140:   return(0);
141: }

143: static PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat)
144: {
146:   PC             pc;
147:   PetscBool      flag;

150:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
151:   KSPGetPC(st->ksp,&pc);
152:   PCGetOperatorsSet(pc,NULL,&flag);
153:   if (flag) {
154:     PCGetOperators(pc,NULL,mat);
155:   } else *mat = NULL;
156:   return(0);
157: }

159: /*@
160:    STPrecondGetMatForPC - Returns the matrix previously set by STPrecondSetMatForPC().

162:    Not Collective, but the Mat is shared by all processors that share the ST

164:    Input Parameter:
165: .  st - the spectral transformation context

167:    Output Parameter:
168: .  mat - the matrix that will be used in constructing the preconditioner or
169:    NULL if no matrix was set by STPrecondSetMatForPC().

171:    Level: advanced

173: .seealso: STPrecondSetMatForPC()
174: @*/
175: PetscErrorCode STPrecondGetMatForPC(ST st,Mat *mat)
176: {

182:   PetscUseMethod(st,"STPrecondGetMatForPC_C",(ST,Mat*),(st,mat));
183:   return(0);
184: }

186: static PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat)
187: {
188:   PC             pc;
189:   Mat            A;
190:   PetscBool      flag;

194:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
195:   KSPGetPC(st->ksp,&pc);
196:   /* Yes, all these lines are needed to safely set mat as the preconditioner
197:      matrix in pc */
198:   PCGetOperatorsSet(pc,&flag,NULL);
199:   if (flag) {
200:     PCGetOperators(pc,&A,NULL);
201:     PetscObjectReference((PetscObject)A);
202:   } else A = NULL;
203:   PetscObjectReference((PetscObject)mat);
204:   PCSetOperators(pc,A,mat);
205:   MatDestroy(&A);
206:   MatDestroy(&mat);
207:   STPrecondSetKSPHasMat(st,PETSC_TRUE);
208:   return(0);
209: }

211: /*@
212:    STPrecondSetMatForPC - Sets the matrix that must be used to build the preconditioner.

214:    Logically Collective on ST and Mat

216:    Input Parameter:
217: +  st - the spectral transformation context
218: -  mat - the matrix that will be used in constructing the preconditioner

220:    Level: advanced

222:    Notes:
223:    This matrix will be passed to the KSP object (via KSPSetOperators) as
224:    the matrix to be used when constructing the preconditioner.
225:    If no matrix is set or mat is set to NULL, A - sigma*B will
226:    be used to build the preconditioner, being sigma the value set by STSetShift().

228: .seealso: STPrecondSetMatForPC(), STSetShift()
229: @*/
230: PetscErrorCode STPrecondSetMatForPC(ST st,Mat mat)
231: {

238:   PetscTryMethod(st,"STPrecondSetMatForPC_C",(ST,Mat),(st,mat));
239:   return(0);
240: }

242: static PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscBool setmat)
243: {
244:   ST_PRECOND *data = (ST_PRECOND*)st->data;

247:   data->setmat = setmat;
248:   return(0);
249: }

251: /*@
252:    STPrecondSetKSPHasMat - Sets a flag indicating that during STSetUp the coefficient
253:    matrix of the KSP linear system (A) must be set to be the same matrix as the
254:    preconditioner (P).

256:    Collective on ST

258:    Input Parameter:
259: +  st - the spectral transformation context
260: -  setmat - the flag

262:    Notes:
263:    In most cases, the preconditioner matrix is used only in the PC object, but
264:    in external solvers this matrix must be provided also as the A-matrix in
265:    the KSP object.

267:    Level: developer

269: .seealso: STPrecondGetKSPHasMat(), STSetShift()
270: @*/
271: PetscErrorCode STPrecondSetKSPHasMat(ST st,PetscBool setmat)
272: {

278:   PetscTryMethod(st,"STPrecondSetKSPHasMat_C",(ST,PetscBool),(st,setmat));
279:   return(0);
280: }

282: static PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscBool *setmat)
283: {
284:   ST_PRECOND *data = (ST_PRECOND*)st->data;

287:   *setmat = data->setmat;
288:   return(0);
289: }

291: /*@
292:    STPrecondGetKSPHasMat - Returns the flag indicating if the coefficient
293:    matrix of the KSP linear system (A) is set to be the same matrix as the
294:    preconditioner (P).

296:    Not Collective

298:    Input Parameter:
299: .  st - the spectral transformation context

301:    Output Parameter:
302: .  setmat - the flag

304:    Level: developer

306: .seealso: STPrecondSetKSPHasMat(), STSetShift()
307: @*/
308: PetscErrorCode STPrecondGetKSPHasMat(ST st,PetscBool *setmat)
309: {

315:   PetscUseMethod(st,"STPrecondGetKSPHasMat_C",(ST,PetscBool*),(st,setmat));
316:   return(0);
317: }

319: PetscErrorCode STDestroy_Precond(ST st)
320: {

324:   PetscFree(st->data);
325:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetMatForPC_C",NULL);
326:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetMatForPC_C",NULL);
327:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",NULL);
328:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",NULL);
329:   return(0);
330: }

332: PETSC_EXTERN PetscErrorCode STCreate_Precond(ST st)
333: {
335:   ST_PRECOND     *ctx;

338:   PetscNewLog(st,&ctx);
339:   st->data = (void*)ctx;

341:   st->ops->getbilinearform = STGetBilinearForm_Default;
342:   st->ops->setup           = STSetUp_Precond;
343:   st->ops->setshift        = STSetShift_Precond;
344:   st->ops->destroy         = STDestroy_Precond;
345:   st->ops->setdefaultksp   = STSetDefaultKSP_Precond;

347:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetMatForPC_C",STPrecondGetMatForPC_Precond);
348:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetMatForPC_C",STPrecondSetMatForPC_Precond);
349:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",STPrecondGetKSPHasMat_Precond);
350:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",STPrecondSetKSPHasMat_Precond);

352:   STPrecondSetKSPHasMat_Precond(st,PETSC_TRUE);
353:   return(0);
354: }