Actual source code: pcksp.c

  1: #define PETSCKSP_DLL

 3:  #include private/pcimpl.h
 4:  #include petscksp.h

  6: typedef struct {
  7:   PetscTruth use_true_matrix;       /* use mat rather than pmat in inner linear solve */
  8:   KSP        ksp;
  9:   PetscInt   its;                   /* total number of iterations KSP uses */
 10: } PC_KSP;

 14: static PetscErrorCode PCApply_KSP(PC pc,Vec x,Vec y)
 15: {
 17:   PetscInt       its;
 18:   PC_KSP         *jac = (PC_KSP*)pc->data;

 21:   KSPSolve(jac->ksp,x,y);
 22:   KSPGetIterationNumber(jac->ksp,&its);
 23:   jac->its += its;
 24:   return(0);
 25: }

 29: static PetscErrorCode PCApplyTranspose_KSP(PC pc,Vec x,Vec y)
 30: {
 32:   PetscInt       its;
 33:   PC_KSP         *jac = (PC_KSP*)pc->data;

 36:   KSPSolveTranspose(jac->ksp,x,y);
 37:   KSPGetIterationNumber(jac->ksp,&its);
 38:   jac->its += its;
 39:   return(0);
 40: }

 44: static PetscErrorCode PCSetUp_KSP(PC pc)
 45: {
 47:   PC_KSP         *jac = (PC_KSP*)pc->data;
 48:   Mat            mat;
 49:   PetscTruth     A;

 52:   KSPSetFromOptions(jac->ksp);
 53:   if (jac->use_true_matrix) mat = pc->mat;
 54:   else                      mat = pc->pmat;

 56:   KSPGetOperatorsSet(jac->ksp,&A,PETSC_NULL);
 57:   if (!A) {
 58:     KSPSetOperators(jac->ksp,mat,pc->pmat,pc->flag);
 59:   }
 60:   KSPSetUp(jac->ksp);
 61:   return(0);
 62: }

 64: /* Default destroy, if it has never been setup */
 67: static PetscErrorCode PCDestroy_KSP(PC pc)
 68: {
 69:   PC_KSP         *jac = (PC_KSP*)pc->data;

 73:   KSPDestroy(jac->ksp);
 74:   PetscFree(jac);
 75:   return(0);
 76: }

 80: static PetscErrorCode PCView_KSP(PC pc,PetscViewer viewer)
 81: {
 82:   PC_KSP         *jac = (PC_KSP*)pc->data;
 84:   PetscTruth     iascii;

 87:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
 88:   if (iascii) {
 89:     if (jac->use_true_matrix) {
 90:       PetscViewerASCIIPrintf(viewer,"Using true matrix (not preconditioner matrix) on inner solve\n");
 91:     }
 92:     PetscViewerASCIIPrintf(viewer,"KSP and PC on KSP preconditioner follow\n");
 93:     PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
 94:   } else {
 95:     SETERRQ1(PETSC_ERR_SUP,"Viewer type %s not supported for this object",((PetscObject)viewer)->type_name);
 96:   }
 97:   PetscViewerASCIIPushTab(viewer);
 98:   KSPView(jac->ksp,viewer);
 99:   PetscViewerASCIIPopTab(viewer);
100:   if (iascii) {
101:     PetscViewerASCIIPrintf(viewer,"---------------------------------\n");
102:   }
103:   return(0);
104: }

108: static PetscErrorCode PCSetFromOptions_KSP(PC pc)
109: {
111:   PetscTruth     flg;

114:   PetscOptionsHead("KSP preconditioner options");
115:     PetscOptionsName("-pc_ksp_true","Use true matrix to define inner linear system, not preconditioner matrix","PCKSPSetUseTrue",&flg);
116:     if (flg) {
117:       PCKSPSetUseTrue(pc);
118:     }
119:   PetscOptionsTail();
120:   return(0);
121: }

123: /* ----------------------------------------------------------------------------------*/

128: PetscErrorCode  PCKSPSetUseTrue_KSP(PC pc)
129: {
130:   PC_KSP   *jac;

133:   jac                  = (PC_KSP*)pc->data;
134:   jac->use_true_matrix = PETSC_TRUE;
135:   return(0);
136: }

142: PetscErrorCode  PCKSPGetKSP_KSP(PC pc,KSP *ksp)
143: {
144:   PC_KSP   *jac;

147:   jac          = (PC_KSP*)pc->data;
148:   *ksp        = jac->ksp;
149:   return(0);
150: }

155: /*@
156:    PCKSPSetUseTrue - Sets a flag to indicate that the true matrix (rather than
157:    the matrix used to define the preconditioner) is used to compute the
158:    residual inside the inner solve.

160:    Collective on PC

162:    Input Parameters:
163: .  pc - the preconditioner context

165:    Options Database Key:
166: .  -pc_ksp_true - Activates PCKSPSetUseTrue()

168:    Note:
169:    For the common case in which the preconditioning and linear 
170:    system matrices are identical, this routine is unnecessary.

172:    Level: advanced

174: .keywords:  PC, KSP, set, true, local, flag

176: .seealso: PCSetOperators(), PCBJacobiSetUseTrueLocal()
177: @*/
178: PetscErrorCode  PCKSPSetUseTrue(PC pc)
179: {
180:   PetscErrorCode ierr,(*f)(PC);

184:   PetscObjectQueryFunction((PetscObject)pc,"PCKSPSetUseTrue_C",(void (**)(void))&f);
185:   if (f) {
186:     (*f)(pc);
187:   }
188:   return(0);
189: }

193: /*@
194:    PCKSPGetKSP - Gets the KSP context for a KSP PC.

196:    Not Collective but KSP returned is parallel if PC was parallel

198:    Input Parameter:
199: .  pc - the preconditioner context

201:    Output Parameters:
202: .  ksp - the PC solver

204:    Notes:
205:    You must call KSPSetUp() before calling PCKSPGetKSP().

207:    Level: advanced

209: .keywords:  PC, KSP, get, context
210: @*/
211: PetscErrorCode  PCKSPGetKSP(PC pc,KSP *ksp)
212: {
213:   PetscErrorCode ierr,(*f)(PC,KSP*);

218:   PetscObjectQueryFunction((PetscObject)pc,"PCKSPGetKSP_C",(void (**)(void))&f);
219:   if (f) {
220:     (*f)(pc,ksp);
221:   }
222:   return(0);
223: }

225: /* ----------------------------------------------------------------------------------*/

227: /*MC
228:      PCKSP -    Defines a preconditioner that can consist of any KSP solver.
229:                  This allows, for example, embedding a Krylov method inside a preconditioner.

231:    Options Database Key:
232: .     -pc_ksp_true - use the matrix that defines the linear system as the matrix for the
233:                     inner solver, otherwise by default it uses the matrix used to construct
234:                     the preconditioner (see PCSetOperators())

236:    Level: intermediate

238:    Concepts: inner iteration

240:    Notes: Using a Krylov method inside another Krylov method can be dangerous (you get divergence or
241:           the incorrect answer) unless you use KSPFGMRES as the other Krylov method


244: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
245:            PCSHELL, PCCOMPOSITE, PCKSPUseTrue(), PCKSPGetKSP()

247: M*/

252: PetscErrorCode  PCCreate_KSP(PC pc)
253: {
255:   const char     *prefix;
256:   PC_KSP         *jac;

259:   PetscNew(PC_KSP,&jac);
260:   PetscLogObjectMemory(pc,sizeof(PC_KSP));
261:   pc->ops->apply              = PCApply_KSP;
262:   pc->ops->applytranspose     = PCApplyTranspose_KSP;
263:   pc->ops->setup              = PCSetUp_KSP;
264:   pc->ops->destroy            = PCDestroy_KSP;
265:   pc->ops->setfromoptions     = PCSetFromOptions_KSP;
266:   pc->ops->view               = PCView_KSP;
267:   pc->ops->applyrichardson    = 0;

269:   pc->data               = (void*)jac;
270:   KSPCreate(pc->comm,&jac->ksp);

272:   PCGetOptionsPrefix(pc,&prefix);
273:   KSPSetOptionsPrefix(jac->ksp,prefix);
274:   KSPAppendOptionsPrefix(jac->ksp,"ksp_");
275:   jac->use_true_matrix = PETSC_FALSE;
276:   jac->its             = 0;

278:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCKSPSetUseTrue_C","PCKSPSetUseTrue_KSP",
279:                     PCKSPSetUseTrue_KSP);
280:   PetscObjectComposeFunctionDynamic((PetscObject)pc,"PCKSPGetKSP_C","PCKSPGetKSP_KSP",
281:                     PCKSPGetKSP_KSP);

283:   return(0);
284: }