Actual source code: precond.c

slepc-3.7.4 2017-05-17
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-2016, 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: static PetscErrorCode STSetDefaultPrecond(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 STSetFromOptions_Precond(PetscOptionItems *PetscOptionsObject,ST st)
 63: {

 67:   STSetDefaultPrecond(st);
 68:   return(0);
 69: }

 73: PetscErrorCode STSetUp_Precond(ST st)
 74: {
 75:   Mat            P;
 76:   PC             pc;
 77:   PetscBool      t0,setmat,destroyP=PETSC_FALSE,builtP;

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

 84:   /* If either pc is none and no matrix has to be set, or pc is shell , exit */
 85:   STSetDefaultPrecond(st);
 86:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
 87:   KSPGetPC(st->ksp,&pc);
 88:   PetscObjectTypeCompare((PetscObject)pc,PCSHELL,&t0);
 89:   if (t0) return(0);
 90:   PetscObjectTypeCompare((PetscObject)pc,PCNONE,&t0);
 91:   STPrecondGetKSPHasMat(st,&setmat);
 92:   if (t0 && !setmat) return(0);

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

 97:   /* If not, create A - shift*B */
 98:   if (P) {
 99:     builtP = PETSC_FALSE;
100:     destroyP = PETSC_TRUE;
101:     PetscObjectReference((PetscObject)P);
102:   } else {
103:     builtP = PETSC_TRUE;

105:     if (!(PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) && st->nmat>1) {
106:       P = st->A[1];
107:       destroyP = PETSC_FALSE;
108:     } else if (st->sigma == 0.0) {
109:       P = st->A[0];
110:       destroyP = PETSC_FALSE;
111:     } else if (PetscAbsScalar(st->sigma) < PETSC_MAX_REAL && st->shift_matrix != ST_MATMODE_SHELL) {
112:       if (st->shift_matrix == ST_MATMODE_INPLACE) {
113:         P = st->A[0];
114:         destroyP = PETSC_FALSE;
115:       } else {
116:         MatDuplicate(st->A[0],MAT_COPY_VALUES,&P);
117:         destroyP = PETSC_TRUE;
118:       }
119:       if (st->nmat>1) {
120:         MatAXPY(P,-st->sigma,st->A[1],st->str);
121:       } else {
122:         MatShift(P,-st->sigma);
123:       }
124:       /* TODO: in case of ST_MATMODE_INPLACE should keep the Hermitian flag of st->A and restore at the end */
125:       STMatSetHermitian(st,P);
126:     } else builtP = PETSC_FALSE;
127:   }

129:   /* If P was not possible to obtain, set pc to PCNONE */
130:   if (!P) {
131:     PCSetType(pc,PCNONE);

133:     /* If some matrix has to be set to ksp, a shell matrix is created */
134:     if (setmat) {
135:       STMatShellCreate(st,-st->sigma,0,NULL,NULL,&P);
136:       STMatSetHermitian(st,P);
137:       destroyP = PETSC_TRUE;
138:     }
139:   }

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

143:   if (destroyP) {
144:     MatDestroy(&P);
145:   } else if (st->shift_matrix == ST_MATMODE_INPLACE && builtP) {
146:     if (st->sigma != 0.0 && PetscAbsScalar(st->sigma) < PETSC_MAX_REAL) {
147:       if (st->nmat>1) {
148:         MatAXPY(st->A[0],st->sigma,st->A[1],st->str);
149:       } else {
150:         MatShift(st->A[0],st->sigma);
151:       }
152:     }
153:   }
154:   return(0);
155: }

159: PetscErrorCode STSetShift_Precond(ST st,PetscScalar newshift)
160: {

164:   /* Nothing to be done if STSetUp has not been called yet */
165:   if (!st->state) return(0);
166:   st->sigma = newshift;
167:   if (st->shift_matrix != ST_MATMODE_SHELL) {
168:     STSetUp_Precond(st);
169:   }
170:   return(0);
171: }

175: static PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat)
176: {
178:   PC             pc;
179:   PetscBool      flag;

182:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
183:   KSPGetPC(st->ksp,&pc);
184:   PCGetOperatorsSet(pc,NULL,&flag);
185:   if (flag) {
186:     PCGetOperators(pc,NULL,mat);
187:   } else *mat = NULL;
188:   return(0);
189: }

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

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

198:    Input Parameter:
199: .  st - the spectral transformation context

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

205:    Level: advanced

207: .seealso: STPrecondSetMatForPC()
208: @*/
209: PetscErrorCode STPrecondGetMatForPC(ST st,Mat *mat)
210: {

216:   PetscUseMethod(st,"STPrecondGetMatForPC_C",(ST,Mat*),(st,mat));
217:   return(0);
218: }

222: static PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat)
223: {
224:   PC             pc;
225:   Mat            A;
226:   PetscBool      flag;

230:   if (!st->ksp) { STGetKSP(st,&st->ksp); }
231:   KSPGetPC(st->ksp,&pc);
232:   /* Yes, all these lines are needed to safely set mat as the preconditioner
233:      matrix in pc */
234:   PCGetOperatorsSet(pc,&flag,NULL);
235:   if (flag) {
236:     PCGetOperators(pc,&A,NULL);
237:     PetscObjectReference((PetscObject)A);
238:   } else A = NULL;
239:   PetscObjectReference((PetscObject)mat);
240:   PCSetOperators(pc,A,mat);
241:   MatDestroy(&A);
242:   MatDestroy(&mat);
243:   STPrecondSetKSPHasMat(st,PETSC_TRUE);
244:   return(0);
245: }

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

252:    Logically Collective on ST and Mat

254:    Input Parameter:
255: +  st - the spectral transformation context
256: -  mat - the matrix that will be used in constructing the preconditioner

258:    Level: advanced

260:    Notes:
261:    This matrix will be passed to the KSP object (via KSPSetOperators) as
262:    the matrix to be used when constructing the preconditioner.
263:    If no matrix is set or mat is set to NULL, A - sigma*B will
264:    be used to build the preconditioner, being sigma the value set by STSetShift().

266: .seealso: STPrecondSetMatForPC(), STSetShift()
267: @*/
268: PetscErrorCode STPrecondSetMatForPC(ST st,Mat mat)
269: {

276:   PetscTryMethod(st,"STPrecondSetMatForPC_C",(ST,Mat),(st,mat));
277:   return(0);
278: }

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

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

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

298:    Collective on ST

300:    Input Parameter:
301: +  st - the spectral transformation context
302: -  setmat - the flag

304:    Notes:
305:    In most cases, the preconditioner matrix is used only in the PC object, but
306:    in external solvers this matrix must be provided also as the A-matrix in
307:    the KSP object.

309:    Level: developer

311: .seealso: STPrecondGetKSPHasMat(), STSetShift()
312: @*/
313: PetscErrorCode STPrecondSetKSPHasMat(ST st,PetscBool setmat)
314: {

320:   PetscTryMethod(st,"STPrecondSetKSPHasMat_C",(ST,PetscBool),(st,setmat));
321:   return(0);
322: }

326: static PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscBool *setmat)
327: {
328:   ST_PRECOND *data = (ST_PRECOND*)st->data;

331:   *setmat = data->setmat;
332:   return(0);
333: }

337: /*@
338:    STPrecondGetKSPHasMat - Returns the flag indicating if the coefficient
339:    matrix of the KSP linear system (A) is set to be the same matrix as the
340:    preconditioner (P).

342:    Not Collective

344:    Input Parameter:
345: .  st - the spectral transformation context

347:    Output Parameter:
348: .  setmat - the flag

350:    Level: developer

352: .seealso: STPrecondSetKSPHasMat(), STSetShift()
353: @*/
354: PetscErrorCode STPrecondGetKSPHasMat(ST st,PetscBool *setmat)
355: {

361:   PetscUseMethod(st,"STPrecondGetKSPHasMat_C",(ST,PetscBool*),(st,setmat));
362:   return(0);
363: }

367: PetscErrorCode STDestroy_Precond(ST st)
368: {

372:   PetscFree(st->data);
373:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetMatForPC_C",NULL);
374:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetMatForPC_C",NULL);
375:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",NULL);
376:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",NULL);
377:   return(0);
378: }

382: PETSC_EXTERN PetscErrorCode STCreate_Precond(ST st)
383: {
385:   ST_PRECOND     *ctx;

388:   PetscNewLog(st,&ctx);
389:   st->data = (void*)ctx;

391:   st->ops->getbilinearform = STGetBilinearForm_Default;
392:   st->ops->setup           = STSetUp_Precond;
393:   st->ops->setshift        = STSetShift_Precond;
394:   st->ops->destroy         = STDestroy_Precond;
395:   st->ops->setfromoptions  = STSetFromOptions_Precond;

397:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetMatForPC_C",STPrecondGetMatForPC_Precond);
398:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetMatForPC_C",STPrecondSetMatForPC_Precond);
399:   PetscObjectComposeFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",STPrecondGetKSPHasMat_Precond);
400:   PetscObjectComposeFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",STPrecondSetKSPHasMat_Precond);

402:   STPrecondSetKSPHasMat_Precond(st,PETSC_TRUE);
403:   return(0);
404: }