Actual source code: precond.c

slepc-3.6.1 2015-09-03
Report Typos and Errors
  1: /*
  2:       Implements the ST class for preconditioned eigenvalue methods.

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

  8:    This file is part of SLEPc.

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

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

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

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

 26: typedef struct {
 27:   PetscBool setmat;
 28: } ST_PRECOND;

 32: PetscErrorCode STSetFromOptions_Precond(PetscOptions *PetscOptionsObject,ST st)
 33: {
 35:   PC             pc;
 36:   PCType         pctype;
 37:   Mat            P;
 38:   PetscBool      t0,t1;
 39:   KSP            ksp;

 42:   STGetKSP(st,&ksp);
 43:   KSPGetPC(ksp,&pc);
 44:   PetscObjectGetType((PetscObject)pc,&pctype);
 45:   STPrecondGetMatForPC(st,&P);
 46:   if (!pctype && st->A && st->A[0]) {
 47:     if (P || st->shift_matrix == ST_MATMODE_SHELL) {
 48:       PCSetType(pc,PCJACOBI);
 49:     } else {
 50:       MatHasOperation(st->A[0],MATOP_DUPLICATE,&t0);
 51:       if (st->nmat>1) {
 52:         MatHasOperation(st->A[0],MATOP_AXPY,&t1);
 53:       } else t1 = PETSC_TRUE;
 54:       PCSetType(pc,(t0 && t1)?PCJACOBI:PCNONE);
 55:     }
 56:   }
 57:   return(0);
 58: }

 62: PetscErrorCode STSetUp_Precond(ST st)
 63: {
 64:   Mat            P;
 65:   PC             pc;
 66:   PetscBool      t0,setmat,destroyP=PETSC_FALSE,builtP;

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

 73:   /* If either pc is none and no matrix has to be set, or pc is shell , exit */
 74:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
 75:   KSPGetPC(st->ksp,&pc);
 76:   PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&t0);
 77:   if (t0) return(0);
 78:   PetscObjectTypeCompare((PetscObject)pc,PCNONE,&t0);
 79:   STPrecondGetKSPHasMat(st,&setmat);
 80:   if (t0 && !setmat) return(0);

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

 85:   /* If not, create A - shift*B */
 86:   if (P) {
 87:     builtP = PETSC_FALSE;
 88:     destroyP = PETSC_TRUE;
 89:     PetscObjectReference((PetscObject)P);
 90:   } else {
 91:     builtP = PETSC_TRUE;

 93:     if (!(PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) && st->nmat>1) {
 94:       P = st->A[1];
 95:       destroyP = PETSC_FALSE;
 96:     } else if (st->sigma == 0.0) {
 97:       P = st->A[0];
 98:       destroyP = PETSC_FALSE;
 99:     } else if (PetscAbsScalar(st->sigma) < PETSC_MAX_REAL && st->shift_matrix != ST_MATMODE_SHELL) {
100:       if (st->shift_matrix == ST_MATMODE_INPLACE) {
101:         P = st->A[0];
102:         destroyP = PETSC_FALSE;
103:       } else {
104:         MatDuplicate(st->A[0],MAT_COPY_VALUES,&P);
105:         destroyP = PETSC_TRUE;
106:       }
107:       if (st->nmat>1) {
108:         MatAXPY(P,-st->sigma,st->A[1],st->str);
109:       } else {
110:         MatShift(P,-st->sigma);
111:       }
112:       /* TODO: in case of ST_MATMODE_INPLACE should keep the Hermitian flag of st->A and restore at the end */
113:       STMatSetHermitian(st,P);
114:     } else builtP = PETSC_FALSE;
115:   }

117:   /* If P was not possible to obtain, set pc to PCNONE */
118:   if (!P) {
119:     PCSetType(pc,PCNONE);

121:     /* If some matrix has to be set to ksp, a shell matrix is created */
122:     if (setmat) {
123:       STMatShellCreate(st,-st->sigma,0,NULL,NULL,&P);
124:       STMatSetHermitian(st,P);
125:       destroyP = PETSC_TRUE;
126:     }
127:   }

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

131:   if (destroyP) {
132:     MatDestroy(&P);
133:   } else if (st->shift_matrix == ST_MATMODE_INPLACE && builtP) {
134:     if (st->sigma != 0.0 && PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) {
135:       if (st->nmat>1) {
136:         MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
137:       } else {
138:         MatShift(st->A[0],st->sigma);
139:       }
140:     }
141:   }
142:   return(0);
143: }

147: PetscErrorCode STSetShift_Precond(ST st,PetscScalar newshift)
148: {

152:   /* Nothing to be done if STSetUp has not been called yet */
153:   if (!st->setupcalled) return(0);
154:   st->sigma = newshift;
155:   if (st->shift_matrix != ST_MATMODE_SHELL) {
156:     STSetUp_Precond(st);
157:   }
158:   return(0);
159: }

163: static PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat)
164: {
166:   PC             pc;
167:   PetscBool      flag;

170:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
171:   KSPGetPC(st->ksp,&pc);
172:   PCGetOperatorsSet(pc,NULL,&flag);
173:   if (flag) {
174:     PCGetOperators(pc,NULL,mat);
175:   } else *mat = NULL;
176:   return(0);
177: }

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

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

186:    Input Parameter:
187: .  st - the spectral transformation context

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

193:    Level: advanced

195: .seealso: STPrecondSetMatForPC()
196: @*/
197: PetscErrorCode STPrecondGetMatForPC(ST st,Mat *mat)
198: {

204:   PetscTryMethod(st,"STPrecondGetMatForPC_C",(ST,Mat*),(st,mat));
205:   return(0);
206: }

210: static PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat)
211: {
212:   PC             pc;
213:   Mat            A;
214:   PetscBool      flag;

218:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
219:   KSPGetPC(st->ksp,&pc);
220:   /* Yes, all these lines are needed to safely set mat as the preconditioner
221:      matrix in pc */
222:   PCGetOperatorsSet(pc,&flag,NULL);
223:   if (flag) {
224:     PCGetOperators(pc,&A,NULL);
225:     PetscObjectReference((PetscObject)A);
226:   } else A = NULL;
227:   PetscObjectReference((PetscObject)mat);
228:   PCSetOperators(pc,A,mat);
229:   MatDestroy(&A);
230:   MatDestroy(&mat);
231:   STPrecondSetKSPHasMat(st,PETSC_TRUE);
232:   return(0);
233: }

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

240:    Logically Collective on ST and Mat

242:    Input Parameter:
243: +  st - the spectral transformation context
244: -  mat - the matrix that will be used in constructing the preconditioner

246:    Level: advanced

248:    Notes:
249:    This matrix will be passed to the KSP object (via KSPSetOperators) as
250:    the matrix to be used when constructing the preconditioner.
251:    If no matrix is set or mat is set to NULL, A - sigma*B will
252:    be used to build the preconditioner, being sigma the value set by STSetShift().

254: .seealso: STPrecondSetMatForPC(), STSetShift()
255: @*/
256: PetscErrorCode STPrecondSetMatForPC(ST st,Mat mat)
257: {

264:   PetscTryMethod(st,"STPrecondSetMatForPC_C",(ST,Mat),(st,mat));
265:   return(0);
266: }

270: static PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscBool setmat)
271: {
272:   ST_PRECOND *data = (ST_PRECOND*)st->data;

275:   data->setmat = setmat;
276:   return(0);
277: }

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

286:    Collective on ST

288:    Input Parameter:
289: +  st - the spectral transformation context
290: -  setmat - the flag

292:    Notes:
293:    In most cases, the preconditioner matrix is used only in the PC object, but
294:    in external solvers this matrix must be provided also as the A-matrix in
295:    the KSP object.

297:    Level: developer

299: .seealso: STPrecondGetKSPHasMat(), STSetShift()
300: @*/
301: PetscErrorCode STPrecondSetKSPHasMat(ST st,PetscBool setmat)
302: {

308:   PetscTryMethod(st,"STPrecondSetKSPHasMat_C",(ST,PetscBool),(st,setmat));
309:   return(0);
310: }

314: static PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscBool *setmat)
315: {
316:   ST_PRECOND *data = (ST_PRECOND*)st->data;

319:   *setmat = data->setmat;
320:   return(0);
321: }

325: /*@
326:    STPrecondGetKSPHasMat - Returns the flag indicating if the coefficient
327:    matrix of the KSP linear system (A) is set to be the same matrix as the
328:    preconditioner (P).

330:    Not Collective

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

335:    Output Parameter:
336: .  setmat - the flag

338:    Level: developer

340: .seealso: STPrecondSetKSPHasMat(), STSetShift()
341: @*/
342: PetscErrorCode STPrecondGetKSPHasMat(ST st,PetscBool *setmat)
343: {

349:   PetscTryMethod(st,"STPrecondGetKSPHasMat_C",(ST,PetscBool*),(st,setmat));
350:   return(0);
351: }

355: PetscErrorCode STDestroy_Precond(ST st)
356: {

360:   PetscFree(st->data);
361:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetMatForPC_C",NULL);
362:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetMatForPC_C",NULL);
363:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",NULL);
364:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",NULL);
365:   return(0);
366: }

370: PETSC_EXTERN PetscErrorCode STCreate_Precond(ST st)
371: {
373:   ST_PRECOND     *ctx;

376:   PetscNewLog(st,&ctx);
377:   st->data = (void*)ctx;

379:   st->ops->getbilinearform = STGetBilinearForm_Default;
380:   st->ops->setup           = STSetUp_Precond;
381:   st->ops->setshift        = STSetShift_Precond;
382:   st->ops->destroy         = STDestroy_Precond;
383:   st->ops->setfromoptions  = STSetFromOptions_Precond;

385:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetMatForPC_C",STPrecondGetMatForPC_Precond);
386:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetMatForPC_C",STPrecondSetMatForPC_Precond);
387:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",STPrecondGetKSPHasMat_Precond);
388:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",STPrecondSetKSPHasMat_Precond);

390:   STPrecondSetKSPHasMat_Precond(st,PETSC_TRUE);
391:   return(0);
392: }