Actual source code: pcksp.c

petsc-3.10.0 2018-09-12
Report Typos and Errors

  2:  #include <petsc/private/pcimpl.h>
  3:  #include <petscksp.h>

  5: typedef struct {
  6:   KSP       ksp;
  7:   PetscInt  its;                    /* total number of iterations KSP uses */
  8: } PC_KSP;


 11: static PetscErrorCode  PCKSPCreateKSP_KSP(PC pc)
 12: {
 14:   const char     *prefix;
 15:   PC_KSP         *jac = (PC_KSP*)pc->data;

 18:   KSPCreate(PetscObjectComm((PetscObject)pc),&jac->ksp);
 19:   KSPSetErrorIfNotConverged(jac->ksp,pc->erroriffailure);
 20:   PetscObjectIncrementTabLevel((PetscObject)jac->ksp,(PetscObject)pc,1);
 21:   PCGetOptionsPrefix(pc,&prefix);
 22:   KSPSetOptionsPrefix(jac->ksp,prefix);
 23:   KSPAppendOptionsPrefix(jac->ksp,"ksp_");
 24:   return(0);
 25: }

 27: static PetscErrorCode PCApply_KSP(PC pc,Vec x,Vec y)
 28: {
 29:   PetscErrorCode     ierr;
 30:   PetscInt           its;
 31:   PC_KSP             *jac = (PC_KSP*)pc->data;
 32:   KSPConvergedReason reason;

 35:   KSPSolve(jac->ksp,x,y);
 36:   KSPGetConvergedReason(jac->ksp,&reason);
 37:   if (reason == KSP_DIVERGED_PCSETUP_FAILED) {
 38:     pc->failedreason = PC_SUBPC_ERROR;
 39:   }
 40:   KSPGetIterationNumber(jac->ksp,&its);
 41:   jac->its += its;
 42:   return(0);
 43: }

 45: static PetscErrorCode PCApplyTranspose_KSP(PC pc,Vec x,Vec y)
 46: {
 48:   PetscInt       its;
 49:   PC_KSP         *jac = (PC_KSP*)pc->data;

 52:   KSPSolveTranspose(jac->ksp,x,y);
 53:   KSPGetIterationNumber(jac->ksp,&its);
 54:   jac->its += its;
 55:   return(0);
 56: }

 58: static PetscErrorCode PCSetUp_KSP(PC pc)
 59: {
 61:   PC_KSP         *jac = (PC_KSP*)pc->data;
 62:   Mat            mat;

 65:   if (!jac->ksp) {
 66:     PCKSPCreateKSP_KSP(pc);
 67:     KSPSetFromOptions(jac->ksp);
 68:   }
 69:   if (pc->useAmat) mat = pc->mat;
 70:   else             mat = pc->pmat;
 71:   KSPSetOperators(jac->ksp,mat,pc->pmat);
 72:   KSPSetUp(jac->ksp);
 73:   return(0);
 74: }

 76: /* Default destroy, if it has never been setup */
 77: static PetscErrorCode PCReset_KSP(PC pc)
 78: {
 79:   PC_KSP         *jac = (PC_KSP*)pc->data;

 83:   KSPReset(jac->ksp);
 84:   return(0);
 85: }

 87: static PetscErrorCode PCDestroy_KSP(PC pc)
 88: {
 89:   PC_KSP         *jac = (PC_KSP*)pc->data;

 93:   KSPDestroy(&jac->ksp);
 94:   PetscObjectComposeFunction((PetscObject)pc,"PCKSPGetKSP_C",NULL);
 95:   PetscObjectComposeFunction((PetscObject)pc,"PCKSPSetKSP_C",NULL);
 96:   PetscFree(pc->data);
 97:   return(0);
 98: }

100: static PetscErrorCode PCView_KSP(PC pc,PetscViewer viewer)
101: {
102:   PC_KSP         *jac = (PC_KSP*)pc->data;
104:   PetscBool      iascii;

107:   if (!jac->ksp) {PCKSPCreateKSP_KSP(pc);}
108:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
109:   if (iascii) {
110:     if (pc->useAmat) {
111:       PetscViewerASCIIPrintf(viewer,"  Using Amat (not Pmat) as operator on inner solve\n");
112:     }
113:     PetscViewerASCIIPrintf(viewer,"  KSP and PC on KSP preconditioner follow\n");
114:     PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n");
115:   }
116:   PetscViewerASCIIPushTab(viewer);
117:   KSPView(jac->ksp,viewer);
118:   PetscViewerASCIIPopTab(viewer);
119:   if (iascii) {
120:     PetscViewerASCIIPrintf(viewer,"  ---------------------------------\n");
121:   }
122:   return(0);
123: }

125: static PetscErrorCode  PCKSPSetKSP_KSP(PC pc,KSP ksp)
126: {
127:   PC_KSP         *jac = (PC_KSP*)pc->data;

131:   PetscObjectReference((PetscObject)ksp);
132:   KSPDestroy(&jac->ksp);
133:   jac->ksp = ksp;
134:   return(0);
135: }

137: /*@
138:    PCKSPSetKSP - Sets the KSP context for a KSP PC.

140:    Collective on PC

142:    Input Parameter:
143: +  pc - the preconditioner context
144: -  ksp - the KSP solver

146:    Notes:
147:    The PC and the KSP must have the same communicator

149:    Level: advanced

151: .keywords:  PC, KSP, set, context
152: @*/
153: PetscErrorCode  PCKSPSetKSP(PC pc,KSP ksp)
154: {

161:   PetscTryMethod(pc,"PCKSPSetKSP_C",(PC,KSP),(pc,ksp));
162:   return(0);
163: }

165: static PetscErrorCode  PCKSPGetKSP_KSP(PC pc,KSP *ksp)
166: {
167:   PC_KSP         *jac = (PC_KSP*)pc->data;

171:   if (!jac->ksp) {PCKSPCreateKSP_KSP(pc);}
172:   *ksp = jac->ksp;
173:   return(0);
174: }

176: /*@
177:    PCKSPGetKSP - Gets the KSP context for a KSP PC.

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

181:    Input Parameter:
182: .  pc - the preconditioner context

184:    Output Parameters:
185: .  ksp - the KSP solver

187:    Notes:
188:    You must call KSPSetUp() before calling PCKSPGetKSP().

190:    If the PC is not a PCKSP object it raises an error

192:    Level: advanced

194: .keywords:  PC, KSP, get, context
195: @*/
196: PetscErrorCode  PCKSPGetKSP(PC pc,KSP *ksp)
197: {

203:   PetscUseMethod(pc,"PCKSPGetKSP_C",(PC,KSP*),(pc,ksp));
204:   return(0);
205: }

207: static PetscErrorCode PCSetFromOptions_KSP(PetscOptionItems *PetscOptionsObject,PC pc)
208: {
209:   PC_KSP         *jac = (PC_KSP*)pc->data;

213:   PetscOptionsHead(PetscOptionsObject,"PC KSP options");
214:   if (jac->ksp) {
215:     KSPSetFromOptions(jac->ksp);
216:    }
217:   PetscOptionsTail();
218:   return(0);
219: }

221: /* ----------------------------------------------------------------------------------*/

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

227:    Options Database Key:
228: .     -pc_use_amat - use the matrix that defines the linear system, Amat as the matrix for the
229:                     inner solver, otherwise by default it uses the matrix used to construct
230:                     the preconditioner, Pmat (see PCSetOperators())

232:    Level: intermediate

234:    Concepts: inner iteration

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

240:    Developer Notes:
241:     If the outer Krylov method has a nonzero initial guess it will compute a new residual based on that initial guess
242:     and pass that as the right hand side into this KSP (and hence this KSP will always have a zero initial guess). For all outer Krylov methods
243:     except Richardson this is neccessary since Krylov methods, even the flexible ones, need to "see" the result of the action of the preconditioner on the
244:     input (current residual) vector, the action of the preconditioner cannot depend also on some other vector (the "initial guess"). For
245:     KSPRICHARDSON it is possible to provide a PCApplyRichardson_PCKSP() that short circuits returning to the KSP object at each iteration to compute the
246:     residual, see for example PCApplyRichardson_SOR(). We do not implement a PCApplyRichardson_PCKSP()  because (1) using a KSP directly inside a Richardson
247:     is not an efficient algorithm anyways and (2) implementing it for its > 1 would essentially require that we implement Richardson (reimplementing the
248:     Richardson code) inside the PCApplyRichardson_PCKSP() leading to duplicate code.

250: .seealso:  PCCreate(), PCSetType(), PCType (for list of available types), PC,
251:            PCSHELL, PCCOMPOSITE, PCSetUseAmat(), PCKSPGetKSP()

253: M*/

255: PETSC_EXTERN PetscErrorCode PCCreate_KSP(PC pc)
256: {
258:   PC_KSP         *jac;

261:   PetscNewLog(pc,&jac);
262:   pc->data = (void*)jac;

264:   PetscMemzero(pc->ops,sizeof(struct _PCOps));
265:   pc->ops->apply           = PCApply_KSP;
266:   pc->ops->applytranspose  = PCApplyTranspose_KSP;
267:   pc->ops->setup           = PCSetUp_KSP;
268:   pc->ops->reset           = PCReset_KSP;
269:   pc->ops->destroy         = PCDestroy_KSP;
270:   pc->ops->setfromoptions  = PCSetFromOptions_KSP;
271:   pc->ops->view            = PCView_KSP;

273:   PetscObjectComposeFunction((PetscObject)pc,"PCKSPGetKSP_C",PCKSPGetKSP_KSP);
274:   PetscObjectComposeFunction((PetscObject)pc,"PCKSPSetKSP_C",PCKSPSetKSP_KSP);
275:   return(0);
276: }