Actual source code: precond.c

  1: /*
  2:       Implements the ST class for preconditioned eigenvalue methods.

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

  8:    This file is part of SLEPc.
  9:       
 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 private/stimpl.h

 26: PetscErrorCode STDestroy_Precond(ST st);
 27: PetscErrorCode STSetFromOptions_Precond(ST st);
 29: PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat);
 30: PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat);
 31: PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscTruth setmat);
 32: PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscTruth *setmat);

 35: typedef struct {
 36:   PetscTruth     setmat;
 37: } ST_PRECOND;


 42: PetscErrorCode SLEPcNotImplemented_Precond(ST st, Vec x, Vec y) {
 43:   SETERRQ(1, "STPrecond does not support some operation. Please, refer to the SLEPc Manual for more information.");
 44: }

 48: PetscErrorCode STSetFromOptions_Precond(ST st)
 49: {
 51:   PC             pc;
 52:   const PCType   pctype;
 53:   Mat            P;
 54:   PetscTruth     t0, t1;


 58:   KSPGetPC(st->ksp, &pc);
 59:   PetscObjectGetType((PetscObject)pc, &pctype);
 60:   STPrecondGetMatForPC(st, &P);
 61:   if (!pctype && st->A) {
 62:     if (P || st->shift_matrix == ST_MATMODE_SHELL) {
 63:       PCSetType(pc, PCJACOBI);
 64:     } else {
 65:       MatHasOperation(st->A, MATOP_DUPLICATE, &t0);
 66:       if (st->B) {
 67:         MatHasOperation(st->A, MATOP_AXPY, &t1);
 68:       } else {
 69:         t1 = PETSC_TRUE;
 70:       }
 71:       PCSetType(pc, (t0 == PETSC_TRUE && t1 == PETSC_TRUE)?
 72:                              PCJACOBI:PCNONE);
 73:     }
 74:   }

 76:   return(0);
 77: }


 82: PetscErrorCode STSetUp_Precond(ST st)
 83: {
 84:   Mat            P;
 85:   PC             pc;
 86:   PetscTruth     t0, setmat, destroyP=PETSC_FALSE, builtP;


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

 94:   /* If pc is none and any matrix has to be set, exit */
 95:   STSetFromOptions_Precond(st);
 96:   KSPGetPC(st->ksp, &pc);
 97:   PetscTypeCompare((PetscObject)pc, PCNONE, &t0);
 98:   STPrecondGetKSPHasMat(st, &setmat);
 99:   if (t0 == PETSC_TRUE && setmat == PETSC_FALSE) return(0);

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

104:   /* If not, create A - shift*B */
105:   if (P) {
106:     builtP = PETSC_FALSE;
107:     destroyP = PETSC_TRUE;
108:   } else {
109:     builtP = PETSC_TRUE;

111:     if (st->shift_matrix == ST_MATMODE_SHELL) {
112:       STMatShellCreate(st,&P);
113:       //TODO: set the apply and apply transpose to st->mat
114:       destroyP = PETSC_TRUE;
115:     } else if (!(PetscAbsScalar(st->sigma) < PETSC_MAX) && st->B) {
116:       P = st->B;
117:       destroyP = PETSC_FALSE;
118:     } else if (st->sigma == 0.0) {
119:       P = st->A;
120:       destroyP = PETSC_FALSE;
121:     } else if (PetscAbsScalar(st->sigma) < PETSC_MAX) {
122:       if (st->shift_matrix == ST_MATMODE_INPLACE) {
123:         P = st->A;
124:         destroyP = PETSC_FALSE;
125:       } else {
126:         MatDuplicate(st->A, MAT_COPY_VALUES, &P);
127:         destroyP = PETSC_TRUE;
128:       }
129:       if (st->B) {
130:         MatAXPY(P, -st->sigma, st->B, st->str);
131:       } else {
132:         MatShift(P, -st->sigma);
133:       }
134:     } else
135:       builtP = PETSC_FALSE;
136:   }

138:   /* If P was not possible to obtain, set pc to PCNONE */
139:   if (!P) {
140:     PCSetType(pc, PCNONE);

142:     /* If some matrix has to be set to ksp, set ksp to KSPPREONLY */
143:     if (setmat == PETSC_TRUE) {
144:       STMatShellCreate(st, &P);
145:       destroyP = PETSC_TRUE;
146:       KSPSetType(st->ksp, KSPPREONLY);
147:     }
148:   }

150:   KSPSetOperators(st->ksp, setmat==PETSC_TRUE?P:PETSC_NULL, P,
151:                          DIFFERENT_NONZERO_PATTERN);

153:   if (destroyP == PETSC_TRUE) {
154:     MatDestroy(P);
155:   } else if (st->shift_matrix == ST_MATMODE_INPLACE && builtP == PETSC_TRUE) {
156:     if (st->sigma != 0.0 && PetscAbsScalar(st->sigma) < PETSC_MAX) {
157:       if (st->B) {
158:         MatAXPY(st->A,st->sigma,st->B,st->str);
159:       } else {
160:         MatShift(st->A,st->sigma);
161:       }
162:     }
163:   }

165:   return(0);
166: }

170: PetscErrorCode STSetShift_Precond(ST st,PetscScalar newshift)
171: {


176:   /* Nothing to be done if STSetUp has not been called yet */
177:   if (!st->setupcalled) return(0);
178: 
179:   st->sigma = newshift;
180:   if (st->shift_matrix != ST_MATMODE_SHELL) {
181:      STSetUp_Precond(st);
182:   }

184:   return(0);
185: }

190: PetscErrorCode STCreate_Precond(ST st)
191: {
193:   ST_PRECOND     *data;


197:   PetscNew(ST_PRECOND, &data);
198:   st->data                 = data;

200:   st->ops->apply           = SLEPcNotImplemented_Precond;
201:   st->ops->getbilinearform = STGetBilinearForm_Default;
202:   st->ops->applytrans      = SLEPcNotImplemented_Precond;
203:   st->ops->postsolve       = PETSC_NULL;
204:   st->ops->backtr          = PETSC_NULL;
205:   st->ops->setup           = STSetUp_Precond;
206:   st->ops->setshift        = STSetShift_Precond;
207:   st->ops->view            = STView_Default;
208:   st->ops->destroy         = STDestroy_Precond;
209:   st->ops->setfromoptions  = STSetFromOptions_Precond;
210: 
211:   st->checknullspace       = PETSC_NULL;

213:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetMatForPC_C","STPrecondGetMatForPC_Precond",STPrecondGetMatForPC_Precond);
214:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetMatForPC_C","STPrecondSetMatForPC_Precond",STPrecondSetMatForPC_Precond);
215:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetKSPHasMat_C","STPrecondGetKSPHasMat_Precond",STPrecondGetKSPHasMat_Precond);
216:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetKSPHasMat_C","STPrecondSetKSPHasMat_Precond",STPrecondSetKSPHasMat_Precond);

218:   STPrecondSetKSPHasMat_Precond(st, PETSC_TRUE);
219:   KSPSetType(st->ksp, KSPPREONLY);

221:   return(0);
222: }

227: PetscErrorCode STDestroy_Precond(ST st)
228: {


233:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetMatForPC_C","",PETSC_NULL);
234:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetMatForPC_C","",PETSC_NULL);
235:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondGetKSPHasMat_C","",PETSC_NULL);
236:   PetscObjectComposeFunctionDynamic((PetscObject)st,"STPrecondSetKSPHasMat_C","",PETSC_NULL);

238:   PetscFree(st->data);

240:   return(0);
241: }

245: /*@
246:    STPrecondGetMatForPC - Gets the matrix previously set by STPrecondSetMatForPC.
247:    This matrix will be passed as parameter in the KSPSetOperator function as
248:    the matrix to be used in constructing the preconditioner.

250:    Collective on ST

252:    Input Parameter:
253: .  st - the spectral transformation context

255:    Output Parameter:
256: .  mat - the matrix that will be used in constructing the preconditioner or
257:    PETSC_NULL if any previous matrix was set by STPrecondSetMatForPC.

259:    Level: advanced

261: .seealso: STPrecondSetMatForPC(), KSPSetOperator()
262: @*/
263: PetscErrorCode STPrecondGetMatForPC(ST st,Mat *mat)
264: {
265:   PetscErrorCode ierr, (*f)(ST,Mat*);

269:   PetscObjectQueryFunction((PetscObject)st,"STPrecondGetMatForPC_C",(void (**)())&f);
270:   if (f) {
271:     (*f)(st,mat);
272:   }
273:   return(0);
274: }

277: PetscErrorCode STPrecondGetMatForPC_Precond(ST st,Mat *mat)
278: {
280:   PC             pc;
281:   PetscTruth     flag;


286:   KSPGetPC(st->ksp, &pc);
287:   PCGetOperatorsSet(pc, PETSC_NULL, &flag);
288:   if (flag == PETSC_TRUE) {
289:     PCGetOperators(pc, PETSC_NULL, mat, PETSC_NULL);
290:   } else
291:     *mat = PETSC_NULL;

293:   return(0);
294: }

299: /*@
300:    STPrecondSetMatForPC - Sets the matrix that will be passed as parameter in
301:    the KSPSetOperators function as the matrix to be used in constructing the
302:    preconditioner. If any matrix is set or mat is PETSC_NULL, A - sigma*B will
303:    be used, being sigma the value set by STSetShift

305:    Collective on ST

307:    Input Parameter:
308: +  st - the spectral transformation context
309: -  mat - the matrix that will be used in constructing the preconditioner

311:    Level: advanced

313: .seealso: STPrecondSetMatForPC(), KSPSetOperators()
314: @*/
315: PetscErrorCode STPrecondSetMatForPC(ST st,Mat mat)
316: {
317:   PetscErrorCode ierr, (*f)(ST,Mat);

322:   PetscObjectQueryFunction((PetscObject)st,"STPrecondGetMatForPC_C",(void (**)())&f);
323:   if (f) {
324:     (*f)(st,mat);
325:   }
326:   return(0);
327: }

330: PetscErrorCode STPrecondSetMatForPC_Precond(ST st,Mat mat)
331: {
332:   PC             pc;


339:   KSPGetPC(st->ksp, &pc);
340:   PCSetOperators(pc, PETSC_NULL, mat, DIFFERENT_NONZERO_PATTERN);

342:   return(0);
343: }


349: /*@
350:    STPrecondSetKSPHasMat - Sets if during the STSetUp the KSP matrix associated
351:    to the linear system is set with the matrix for building the preconditioner.

353:    Collective on ST

355:    Input Parameter:
356: +  st - the spectral transformation context
357: -  setmat - if true, the KSP matrix associated to linear system is set with
358:    the matrix for building the preconditioner

360:    Level: developer

362: .seealso: STPrecondGetKSPHasMat(), TSetShift(), KSPSetOperators()
363: @*/
364: PetscErrorCode STPrecondSetKSPHasMat(ST st,PetscTruth setmat)
365: {
366:   PetscErrorCode ierr, (*f)(ST,PetscTruth);

370:   PetscObjectQueryFunction((PetscObject)st,"STPrecondSetKSPHasMat_C",(void (**)())&f);
371:   if (f) {
372:     (*f)(st,setmat);
373:   }
374:   return(0);
375: }

379: /*@
380:    STPrecondGetKSPHasMat - Gets if during the STSetUp the KSP matrix associated
381:    to the linear system is set with the matrix for building the preconditioner.

383:    Collective on ST

385:    Input Parameter:
386: .  st - the spectral transformation context

388:    Output Parameter:
389: .  setmat - if true, the KSP matrix associated to linear system is set with
390:    the matrix for building the preconditioner

392:    Level: developer

394: .seealso: STPrecondSetKSPHasMat(), STSetShift(), KSPSetOperators()
395: @*/
396: PetscErrorCode STPrecondGetKSPHasMat(ST st,PetscTruth *setmat)
397: {
398:   PetscErrorCode ierr, (*f)(ST,PetscTruth*);

402:   PetscObjectQueryFunction((PetscObject)st,"STPrecondGetKSPHasMat_C",(void (**)())&f);
403:   if (f) {
404:     (*f)(st,setmat);
405:   }
406:   return(0);
407: }

410: PetscErrorCode STPrecondSetKSPHasMat_Precond(ST st,PetscTruth setmat)
411: {
412:   ST_PRECOND     *data = (ST_PRECOND*)st->data;



418:   data->setmat = setmat;

420:   return(0);
421: }

425: PetscErrorCode STPrecondGetKSPHasMat_Precond(ST st,PetscTruth *setmat)
426: {
427:   ST_PRECOND     *data = (ST_PRECOND*)st->data;



433:   *setmat = data->setmat;

435:   return(0);
436: }