Actual source code: gcr.c

petsc-3.11.4 2019-09-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:       KSPCheckNorm(ksp,norm_r);
 65:     }
 66:     /* update the local counter and the global counter */
 67:     ksp->its++;
 68:     ksp->rnorm = norm_r;

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

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

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

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

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

101:   /* compute initial residual */
102:   KSP_MatMult(ksp,A, x, r);
103:   VecAYPX(r, -1.0, b); /* r = b - A x  */
104:   VecNorm(r, NORM_2, &norm_r);
105:   KSPCheckNorm(ksp,norm_r);
106:   ksp->its    = 0;
107:   ksp->rnorm0 = norm_r;

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

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

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

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

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


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

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

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

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

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

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

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

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

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

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

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

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

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

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

219:  Logically Collective on KSP

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

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

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

235:  Level: intermediate

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

240:  .seealso: KSPGCRModifyPCNoChange()

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

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

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

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

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

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

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

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

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

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

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

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

309:    Level: beginner

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

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

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

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

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

331:     Contributed by Dave May

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


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

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

348:   PetscNewLog(ksp,&ctx);

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

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

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

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