Actual source code: pcksp.c

petsc-3.8.3 2017-12-09
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) {PCKSPCreateKSP_KSP(pc);}
 66:   KSPSetFromOptions(jac->ksp);
 67:   if (pc->useAmat) mat = pc->mat;
 68:   else             mat = pc->pmat;

 70:   KSPSetOperators(jac->ksp,mat,pc->pmat);
 71:   KSPSetUp(jac->ksp);
 72:   return(0);
 73: }

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

 82:   if (jac->ksp) {KSPReset(jac->ksp);}
 83:   return(0);
 84: }

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

 92:   PCReset_KSP(pc);
 93:   KSPDestroy(&jac->ksp);
 94:   PetscFree(pc->data);
 95:   return(0);
 96: }

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

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

123: static PetscErrorCode  PCKSPGetKSP_KSP(PC pc,KSP *ksp)
124: {
125:   PC_KSP         *jac = (PC_KSP*)pc->data;

129:   if (!jac->ksp) {PCKSPCreateKSP_KSP(pc);}
130:   *ksp = jac->ksp;
131:   return(0);
132: }

134: /*@
135:    PCKSPGetKSP - Gets the KSP context for a KSP PC.

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

139:    Input Parameter:
140: .  pc - the preconditioner context

142:    Output Parameters:
143: .  ksp - the PC solver

145:    Notes:
146:    You must call KSPSetUp() before calling PCKSPGetKSP().

148:    If the PC is not a PCKSP object then a NULL is returned

150:    Level: advanced

152: .keywords:  PC, KSP, get, context
153: @*/
154: PetscErrorCode  PCKSPGetKSP(PC pc,KSP *ksp)
155: {

161:   *ksp = NULL;
162:   PetscTryMethod(pc,"PCKSPGetKSP_C",(PC,KSP*),(pc,ksp));
163:   return(0);
164: }

166: static PetscErrorCode PCSetFromOptions_KSP(PetscOptionItems *PetscOptionsObject,PC pc)
167: {
168:   PC_KSP         *jac = (PC_KSP*)pc->data;

172:   PetscOptionsHead(PetscOptionsObject,"PC KSP options");
173:   if (jac->ksp) {
174:     KSPSetFromOptions(jac->ksp);
175:    }
176:   PetscOptionsTail();
177:   return(0);
178: }

180: /* ----------------------------------------------------------------------------------*/

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

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

191:    Level: intermediate

193:    Concepts: inner iteration

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

198:    Developer Notes: If the outer Krylov method has a nonzero initial guess it will compute a new residual based on that initial guess
199:     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
200:     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
201:     input (current residual) vector, the action of the preconditioner cannot depend also on some other vector (the "initial guess"). For 
202:     KSPRICHARDSON it is possible to provide a PCApplyRichardson_PCKSP() that short circuits returning to the KSP object at each iteration to compute the 
203:     residual, see for example PCApplyRichardson_SOR(). We do not implement a PCApplyRichardson_PCKSP()  because (1) using a KSP directly inside a Richardson 
204:     is not an efficient algorithm anyways and (2) implementing it for its > 1 would essentially require that we implement Richardson (reimplementing the 
205:     Richardson code) inside the PCApplyRichardson_PCKSP() leading to duplicate code.

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

210: M*/

212: PETSC_EXTERN PetscErrorCode PCCreate_KSP(PC pc)
213: {
215:   PC_KSP         *jac;

218:   PetscNewLog(pc,&jac);

220:   pc->ops->apply           = PCApply_KSP;
221:   pc->ops->applytranspose  = PCApplyTranspose_KSP;
222:   pc->ops->setup           = PCSetUp_KSP;
223:   pc->ops->reset           = PCReset_KSP;
224:   pc->ops->destroy         = PCDestroy_KSP;
225:   pc->ops->setfromoptions  = PCSetFromOptions_KSP;
226:   pc->ops->view            = PCView_KSP;
227:   pc->ops->applyrichardson = 0;

229:   pc->data = (void*)jac;


232:   jac->its             = 0;
233:   PetscObjectComposeFunction((PetscObject)pc,"PCKSPGetKSP_C",PCKSPGetKSP_KSP);
234:   return(0);
235: }