Actual source code: gcr.c

petsc-3.10.5 2019-03-28
Report Typos and Errors

  2:  #include <petsc/private/kspimpl.h>

  4: typedef struct {
  5:   PetscInt    restart;
  6:   PetscInt    n_restarts;
  7:   PetscScalar *val;
  8:   Vec         *VV, *SS;
  9:   Vec         R;

 11:   PetscErrorCode (*modifypc)(KSP,PetscInt,PetscReal,void*);  /* function to modify the preconditioner*/
 12:   PetscErrorCode (*modifypc_destroy)(void*);                 /* function to destroy the user context for the modifypc function */

 14:   void *modifypc_ctx;                                        /* user defined data for the modifypc function */
 15: } KSP_GCR;

 17: static PetscErrorCode KSPSolve_GCR_cycle(KSP ksp)
 18: {
 19:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
 21:   PetscScalar    r_dot_v;
 22:   Mat            A, B;
 23:   PC             pc;
 24:   Vec            s,v,r;
 25:   /*
 26:      The residual norm will not be computed when ksp->its > ksp->chknorm hence need to initialize norm_r with some dummy value
 27:   */
 28:   PetscReal      norm_r = 0.0,nrm;
 29:   PetscInt       k, i, restart;
 30:   Vec            x;

 33:   restart = ctx->restart;
 34:   KSPGetPC(ksp, &pc);
 35:   KSPGetOperators(ksp, &A, &B);

 37:   x = ksp->vec_sol;
 38:   r = ctx->R;

 40:   for (k=0; k<restart; k++) {
 41:     v = ctx->VV[k];
 42:     s = ctx->SS[k];
 43:     if (ctx->modifypc) {
 44:       (*ctx->modifypc)(ksp,ksp->its,ksp->rnorm,ctx->modifypc_ctx);
 45:     }

 47:     KSP_PCApply(ksp, r, s); /* s = B^{-1} r */
 48:     KSP_MatMult(ksp,A, s, v);  /* v = A s */

 50:     VecMDot(v,k, ctx->VV, ctx->val);
 51:     for (i=0; i<k; i++) ctx->val[i] = -ctx->val[i];
 52:     VecMAXPY(v,k,ctx->val,ctx->VV); /* v = v - sum_{i=0}^{k-1} alpha_i v_i */
 53:     VecMAXPY(s,k,ctx->val,ctx->SS); /* s = s - sum_{i=0}^{k-1} alpha_i s_i */

 55:     VecDotNorm2(r,v,&r_dot_v,&nrm);
 56:     nrm     = PetscSqrtReal(nrm);
 57:     r_dot_v = r_dot_v/nrm;
 58:     VecScale(v, 1.0/nrm);
 59:     VecScale(s, 1.0/nrm);
 60:     VecAXPY(x,  r_dot_v, s);
 61:     VecAXPY(r, -r_dot_v, v);
 62:     if (ksp->its > ksp->chknorm) {
 63:       VecNorm(r, NORM_2, &norm_r);
 64:     }
 65:     /* update the local counter and the global counter */
 66:     ksp->its++;
 67:     ksp->rnorm = norm_r;

 69:     KSPLogResidualHistory(ksp,norm_r);
 70:     KSPMonitor(ksp,ksp->its,norm_r);

 72:     if (ksp->its-1 > ksp->chknorm) {
 73:       (*ksp->converged)(ksp,ksp->its,norm_r,&ksp->reason,ksp->cnvP);
 74:       if (ksp->reason) break;
 75:     }

 77:     if (ksp->its >= ksp->max_it) {
 78:       ksp->reason = KSP_CONVERGED_ITS;
 79:       break;
 80:     }
 81:   }
 82:   ctx->n_restarts++;
 83:   return(0);
 84: }

 86: static PetscErrorCode KSPSolve_GCR(KSP ksp)
 87: {
 88:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
 90:   Mat            A, B;
 91:   Vec            r,b,x;
 92:   PetscReal      norm_r;

 95:   KSPGetOperators(ksp, &A, &B);
 96:   x    = ksp->vec_sol;
 97:   b    = ksp->vec_rhs;
 98:   r    = ctx->R;

100:   /* compute initial residual */
101:   KSP_MatMult(ksp,A, x, r);
102:   VecAYPX(r, -1.0, b); /* r = b - A x  */
103:   VecNorm(r, NORM_2, &norm_r);

105:   ksp->its    = 0;
106:   ksp->rnorm0 = norm_r;

108:   KSPLogResidualHistory(ksp,ksp->rnorm0);
109:   KSPMonitor(ksp,ksp->its,ksp->rnorm0);
110:   (*ksp->converged)(ksp,ksp->its,ksp->rnorm0,&ksp->reason,ksp->cnvP);
111:   if (ksp->reason) return(0);

113:   do {
114:     KSPSolve_GCR_cycle(ksp);
115:     if (ksp->reason) break; /* catch case when convergence occurs inside the cycle */
116:   } while (ksp->its < ksp->max_it);

118:   if (ksp->its >= ksp->max_it) ksp->reason = KSP_DIVERGED_ITS;
119:   return(0);
120: }

122: static PetscErrorCode KSPView_GCR(KSP ksp, PetscViewer viewer)
123: {
124:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
126:   PetscBool      iascii;

129:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
130:   if (iascii) {
131:     PetscViewerASCIIPrintf(viewer,"  restart = %D \n", ctx->restart);
132:     PetscViewerASCIIPrintf(viewer,"  restarts performed = %D \n", ctx->n_restarts);
133:   }
134:   return(0);
135: }


138: static PetscErrorCode KSPSetUp_GCR(KSP ksp)
139: {
140:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
142:   Mat            A;
143:   PetscBool      diagonalscale;

146:   PCGetDiagonalScale(ksp->pc,&diagonalscale);
147:   if (diagonalscale) SETERRQ1(PetscObjectComm((PetscObject)ksp),PETSC_ERR_SUP,"Krylov method %s does not support diagonal scaling",((PetscObject)ksp)->type_name);

149:   KSPGetOperators(ksp, &A, NULL);
150:   MatCreateVecs(A, &ctx->R, NULL);
151:   VecDuplicateVecs(ctx->R, ctx->restart, &ctx->VV);
152:   VecDuplicateVecs(ctx->R, ctx->restart, &ctx->SS);

154:   PetscMalloc1(ctx->restart, &ctx->val);
155:   return(0);
156: }

158: static PetscErrorCode KSPReset_GCR(KSP ksp)
159: {
161:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;

164:   VecDestroy(&ctx->R);
165:   VecDestroyVecs(ctx->restart,&ctx->VV);
166:   VecDestroyVecs(ctx->restart,&ctx->SS);
167:   if (ctx->modifypc_destroy) {
168:     (*ctx->modifypc_destroy)(ctx->modifypc_ctx);
169:   }
170:   PetscFree(ctx->val);
171:   return(0);
172: }

174: static PetscErrorCode KSPDestroy_GCR(KSP ksp)
175: {

179:   KSPReset_GCR(ksp);
180:   KSPDestroyDefault(ksp);
181:   return(0);
182: }

184: static PetscErrorCode KSPSetFromOptions_GCR(PetscOptionItems *PetscOptionsObject,KSP ksp)
185: {
187:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
188:   PetscInt       restart;
189:   PetscBool      flg;

192:   PetscOptionsHead(PetscOptionsObject,"KSP GCR options");
193:   PetscOptionsInt("-ksp_gcr_restart","Number of Krylov search directions","KSPGCRSetRestart",ctx->restart,&restart,&flg);
194:   if (flg) { KSPGCRSetRestart(ksp,restart); }
195:   PetscOptionsTail();
196:   return(0);
197: }

200: typedef PetscErrorCode (*KSPGCRModifyPCFunction)(KSP,PetscInt,PetscReal,void*);
201: typedef PetscErrorCode (*KSPGCRDestroyFunction)(void*);

203: static PetscErrorCode  KSPGCRSetModifyPC_GCR(KSP ksp,KSPGCRModifyPCFunction function,void *data,KSPGCRDestroyFunction destroy)
204: {
205:   KSP_GCR *ctx = (KSP_GCR*)ksp->data;

209:   ctx->modifypc         = function;
210:   ctx->modifypc_destroy = destroy;
211:   ctx->modifypc_ctx     = data;
212:   return(0);
213: }

215: /*@C
216:  KSPGCRSetModifyPC - Sets the routine used by GCR to modify the preconditioner.

218:  Logically Collective on KSP

220:  Input Parameters:
221:  +  ksp      - iterative context obtained from KSPCreate()
222:  .  function - user defined function to modify the preconditioner
223:  .  ctx      - user provided contex for the modify preconditioner function
224:  -  destroy  - the function to use to destroy the user provided application context.

226:  Calling Sequence of function:
227:   PetscErrorCode function (KSP ksp, PetscInt n, PetscReal rnorm, void *ctx)

229:  ksp   - iterative context
230:  n     - the total number of GCR iterations that have occurred
231:  rnorm - 2-norm residual value
232:  ctx   - the user provided application context

234:  Level: intermediate

236:  Notes:
237:  The default modifypc routine is KSPGCRModifyPCNoChange()

239:  .seealso: KSPGCRModifyPCNoChange()

241:  @*/
242: PetscErrorCode  KSPGCRSetModifyPC(KSP ksp,PetscErrorCode (*function)(KSP,PetscInt,PetscReal,void*),void *data,PetscErrorCode (*destroy)(void*))
243: {

247:   PetscUseMethod(ksp,"KSPGCRSetModifyPC_C",(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *data,PetscErrorCode (*)(void*)),(ksp,function,data,destroy));
248:   return(0);
249: }

251: static PetscErrorCode KSPGCRSetRestart_GCR(KSP ksp,PetscInt restart)
252: {
253:   KSP_GCR *ctx;

256:   ctx          = (KSP_GCR*)ksp->data;
257:   ctx->restart = restart;
258:   return(0);
259: }

261: PetscErrorCode  KSPGCRSetRestart(KSP ksp, PetscInt restart)
262: {

266:   PetscTryMethod(ksp,"KSPGCRSetRestart_C",(KSP,PetscInt),(ksp,restart));
267:   return(0);
268: }

270: static PetscErrorCode  KSPBuildSolution_GCR(KSP ksp, Vec v, Vec *V)
271: {
273:   Vec            x;

276:   x = ksp->vec_sol;
277:   if (v) {
278:     VecCopy(x, v);
279:     if (V) *V = v;
280:   } else if (V) {
281:     *V = ksp->vec_sol;
282:   }
283:   return(0);
284: }

286: static PetscErrorCode  KSPBuildResidual_GCR(KSP ksp, Vec t, Vec v, Vec *V)
287: {
289:   KSP_GCR        *ctx;

292:   ctx = (KSP_GCR*)ksp->data;
293:   if (v) {
294:     VecCopy(ctx->R, v);
295:     if (V) *V = v;
296:   } else if (V) {
297:     *V = ctx->R;
298:   }
299:   return(0);
300: }

302: /*MC
303:      KSPGCR - Implements the preconditioned Generalized Conjugate Residual method.

305:    Options Database Keys:
306: .   -ksp_gcr_restart <restart> - the number of stored vectors to orthogonalize against

308:    Level: beginner

310:     Notes:
311:     The GCR Krylov method supports non-symmetric matrices and permits the use of a preconditioner
312:            which may vary from one iteration to the next. Users can can define a method to vary the
313:            preconditioner between iterates via KSPGCRSetModifyPC().

315:            Restarts are solves with x0 not equal to zero. When a restart occurs, the initial starting
316:            solution is given by the current estimate for x which was obtained by the last restart
317:            iterations of the GCR algorithm.

319:            Unlike GMRES and FGMRES, when using GCR, the solution and residual vector can be directly accessed at any iterate,
320:            with zero computational cost, via a call to KSPBuildSolution() and KSPBuildResidual() respectively.

322:            This implementation of GCR will only apply the stopping condition test whenever ksp->its > ksp->chknorm,
323:            where ksp->chknorm is specified via the command line argument -ksp_check_norm_iteration or via
324:            the function KSPSetCheckNormIteration(). Hence the residual norm reported by the monitor and stored
325:            in the residual history will be listed as 0.0 before this iteration. It is actually not 0.0; just not calculated.

327:            The method implemented requires the storage of 2 x restart + 1 vectors, twice as much as GMRES.
328:            Support only for right preconditioning.

330:     Contributed by Dave May

332:     References:
333: .          1. - S. C. Eisenstat, H. C. Elman, and H. C. Schultz. Variational iterative methods for
334:            nonsymmetric systems of linear equations. SIAM J. Numer. Anal., 20, 1983


337: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
338:            KSPGCRSetRestart(), KSPGCRSetModifyPC(), KSPGMRES, KSPFGMRES

340: M*/
341: PETSC_EXTERN PetscErrorCode KSPCreate_GCR(KSP ksp)
342: {
344:   KSP_GCR        *ctx;

347:   PetscNewLog(ksp,&ctx);

349:   ctx->restart    = 30;
350:   ctx->n_restarts = 0;
351:   ksp->data       = (void*)ctx;

353:   KSPSetSupportedNorm(ksp,KSP_NORM_UNPRECONDITIONED,PC_RIGHT,3);

355:   ksp->ops->setup          = KSPSetUp_GCR;
356:   ksp->ops->solve          = KSPSolve_GCR;
357:   ksp->ops->reset          = KSPReset_GCR;
358:   ksp->ops->destroy        = KSPDestroy_GCR;
359:   ksp->ops->view           = KSPView_GCR;
360:   ksp->ops->setfromoptions = KSPSetFromOptions_GCR;
361:   ksp->ops->buildsolution  = KSPBuildSolution_GCR;
362:   ksp->ops->buildresidual  = KSPBuildResidual_GCR;

364:   PetscObjectComposeFunction((PetscObject)ksp, "KSPGCRSetRestart_C",KSPGCRSetRestart_GCR);
365:   PetscObjectComposeFunction((PetscObject)ksp,"KSPGCRSetModifyPC_C",KSPGCRSetModifyPC_GCR);
366:   return(0);
367: }