Actual source code: gcr.c

 2:  #include petscksp.h
 3:  #include private/kspimpl.h

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

 12:   PetscErrorCode (*modifypc)(KSP,PetscInt,PetscReal,void*);  /* function to modify the preconditioner*/
 13:   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;

 19: PetscErrorCode KSPSolve_GCR_cycle( KSP ksp )
 20: {
 21:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
 23:   PetscScalar    nrm,r_dot_v;
 24:   Mat            A, B;
 25:   PC             pc;
 26:   Vec            s,v,r;
 27:   PetscReal      norm_r;
 28:   PetscInt       k, i, restart;
 29:   Vec            b,x;
 30:   PetscReal      res;
 31: 
 33:   restart = ctx->restart;
 34:   KSPGetPC( ksp, &pc );
 35:   KSPGetOperators( ksp, &A, &B, 0 );
 36: 
 37:   x = ksp->vec_sol;
 38:   b = ksp->vec_rhs;
 39:   r = ctx->R;

 41:   for ( k=0; k<restart; k++ ) {
 42:     v = ctx->VV[k];
 43:     s = ctx->SS[k];
 44:     if (ctx->modifypc) {
 45:       (*ctx->modifypc)(ksp,ksp->its,ksp->rnorm,ctx->modifypc_ctx);
 46:     }
 47: 
 48:     PCApply( pc, r, s );  /* s = B^{-1} r */
 49:     MatMult( A, s, v );   /* v = A s */
 50: 
 51:     VecMDot( v,k, ctx->VV, ctx->val );
 52: 
 53:     for( i=0; i<=k; i++ ) {
 54:       VecAXPY( v, -ctx->val[i], ctx->VV[i] );  /* v = v - alpha_i v_i */
 55:       VecAXPY( s, -ctx->val[i], ctx->SS[i] );  /* s = s - alpha_i s_i */
 56:     }
 57: 
 58:     VecDotNorm2(r,v,&r_dot_v,&nrm);
 59:     nrm     = sqrt(nrm);
 60:     r_dot_v = r_dot_v/nrm;
 61:     VecScale( v, 1.0/nrm );
 62:     VecScale( s, 1.0/nrm );
 63:     VecAXPY( x,  r_dot_v, s );
 64:     VecAXPY( r, -r_dot_v, v );
 65:     if (ksp->its > ksp->chknorm  ) {
 66:       VecNorm( r, NORM_2, &norm_r );
 67:     }
 68:     /* update the local counter and the global counter */
 69:     ksp->its++;
 70:     res = norm_r;
 71:     ksp->rnorm = res;
 72: 
 73:     KSPLogResidualHistory(ksp,res);
 74:     KSPMonitor(ksp,ksp->its,res);
 75: 
 76:     if( ksp->its > ksp->chknorm  ) {
 77:       (*ksp->converged)(ksp,ksp->its,res,&ksp->reason,ksp->cnvP);
 78:       if (ksp->reason) break;
 79:     }
 80: 
 81:     if( ksp->its >= ksp->max_it ) {
 82:       ksp->reason = KSP_CONVERGED_ITS;
 83:       break;
 84:     }
 85:   }
 86:   ctx->n_restarts++;
 87:   return(0);
 88: }

 92: PetscErrorCode KSPSolve_GCR( KSP ksp )
 93: {
 94:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
 96:   Mat            A, B;
 97:   Vec            r,b,x;
 98:   PetscReal      norm_r;
 99: 
101:   KSPGetOperators( ksp, &A, &B, PETSC_NULL );
102:   x = ksp->vec_sol;
103:   b = ksp->vec_rhs;
104:   r = ctx->R;
105: 
106:   /* compute initial residual */
107:   MatMult( A, x, r );
108:   VecAYPX( r, -1.0, b ); /* r = b - A x  */
109:   VecNorm( r, NORM_2, &norm_r );
110: 
111:   ksp->its = 0;
112:   ksp->rnorm0 = norm_r;
113: 
114:   KSPLogResidualHistory(ksp,ksp->rnorm0);
115:   KSPMonitor(ksp,ksp->its,ksp->rnorm0);
116:   (*ksp->converged)(ksp,ksp->its,ksp->rnorm0,&ksp->reason,ksp->cnvP);
117:   if (ksp->reason) return(0);
118: 
119:   do {
120:     KSPSolve_GCR_cycle( ksp );
121:     if (ksp->reason) break; /* catch case when convergence occurs inside the cycle */
122:   } while( ksp->its < ksp->max_it );
123:   if( ksp->its >= ksp->max_it) {
124:     ksp->reason = KSP_DIVERGED_ITS;
125:   }
126:   return(0);
127: }

131: PetscErrorCode KSPView_GCR( KSP ksp, PetscViewer viewer )
132: {
133:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
135:   PetscTruth     iascii;

138:   PetscTypeCompare((PetscObject)viewer,PETSC_VIEWER_ASCII,&iascii);
139:   if (iascii) {
140:     PetscViewerASCIIPrintf(viewer,"  GCR: restart = %D \n", ctx->restart );
141:     PetscViewerASCIIPrintf(viewer,"  GCR: restarts performed = %D \n", ctx->n_restarts );
142:   }
143:   return(0);
144: }


149: PetscErrorCode KSPSetUp_GCR( KSP ksp )
150: {
151:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
153:   Mat            A;
154: 
156:   if (ksp->pc_side == PC_LEFT) {SETERRQ(PETSC_ERR_SUP,"No left preconditioning for GCR");}
157:   else if (ksp->pc_side == PC_SYMMETRIC) {SETERRQ(PETSC_ERR_SUP,"No symmetric preconditioning for GCR");}

159:   KSPGetOperators( ksp, &A, 0, 0 );
160:   MatGetVecs( A, &ctx->R, PETSC_NULL );
161:   VecDuplicateVecs( ctx->R, ctx->restart, &ctx->VV );
162:   VecDuplicateVecs( ctx->R, ctx->restart, &ctx->SS );
163: 
164:   PetscMalloc( sizeof(PetscScalar)*ctx->restart, &ctx->val );
165:   return(0);
166: }

170: PetscErrorCode KSPDestroy_GCR( KSP ksp )
171: {
173:   KSP_GCR        *ctx = (KSP_GCR*)ksp->data;
174: 
176:   VecDestroy( ctx->R );
177:   VecDestroyVecs( ctx->VV, ctx->restart );
178:   VecDestroyVecs( ctx->SS, ctx->restart );
179:   PetscFree( ctx->val );
180:   if (ctx->modifypc_destroy) {
181:     (*ctx->modifypc_destroy)(ctx->modifypc_ctx);
182:   }
183:   KSPDefaultDestroy(ksp);
184:   return(0);
185: }

189: PetscErrorCode KSPSetFromOptions_GCR(KSP ksp)
190: {
191:   PetscErrorCode  ierr;
192:   KSP_GCR         *ctx = (KSP_GCR *)ksp->data;
193:   PetscInt        restart;
194:   PetscTruth      flg;
195: 
197:   PetscOptionsHead("KSP GCR options");
198:     PetscOptionsInt("-ksp_gcr_restart","Number of Krylov search directions","KSPGCRSetRestart",ctx->restart,&restart,&flg);
199:     if (flg) { KSPGCRSetRestart(ksp,restart); }
200:   PetscOptionsTail();
201:   return(0);
202: }

207: PetscErrorCode  KSPGCRSetModifyPC_GCR(KSP ksp,PetscErrorCode (*function)(KSP,PetscInt,PetscReal,void*),void *data,PetscErrorCode (*destroy)(void*))
208: {
209:   KSP_GCR         *ctx = (KSP_GCR *)ksp->data;
210: 
213:   ctx->modifypc         = function;
214:   ctx->modifypc_destroy = destroy;
215:   ctx->modifypc_ctx     = data;
216:   return(0);
217: }

222: /*@C
223:  KSPGCRSetModifyPC - Sets the routine used by GCR to modify the preconditioner.
224:  
225:  Collective on KSP
226:  
227:  Input Parameters:
228:  +  ksp      - iterative context obtained from KSPCreate()
229:  .  function - user defined function to modify the preconditioner
230:  .  ctx      - user provided contex for the modify preconditioner function
231:  -  destroy  - the function to use to destroy the user provided application context.
232:  
233:  Calling Sequence of function:
234:   PetscErrorCode function ( KSP ksp, PetscInt n, PetscReal rnorm, void *ctx )
235:  
236:  ksp   - iterative context 
237:  n     - the total number of GCR iterations that have occurred     
238:  rnorm - 2-norm residual value
239:  ctx   - the user provided application context
240:  
241:  Level: intermediate
242:  
243:  Notes:
244:  The default modifypc routine is KSPGCRModifyPCNoChange()
245:  
246:  .seealso: KSPGCRModifyPCNoChange()
247:  
248:  @*/
249: PetscErrorCode  KSPGCRSetModifyPC(KSP ksp,PetscErrorCode (*function)(KSP,PetscInt,PetscReal,void*),void *data,PetscErrorCode (*destroy)(void*))
250: {
251:   PetscErrorCode ierr,(*f)(KSP,PetscErrorCode (*)(KSP,PetscInt,PetscReal,void*),void *data,PetscErrorCode (*)(void*));
252: 
254:   PetscObjectQueryFunction((PetscObject)ksp,"KSPGCRSetModifyPC_C",(void (**)(void))&f);
255:   if (f) {
256:     (*f)(ksp,function,data,destroy);
257:   }
258:   return(0);
259: }

263: PetscErrorCode KSPGCRSetRestart_GCR(KSP ksp,PetscInt restart)
264: {
265:   KSP_GCR *ctx;
266: 
268:   ctx = (KSP_GCR *)ksp->data;
269:   ctx->restart = restart;
270:   return(0);
271: }

275: PetscErrorCode  KSPGCRSetRestart(KSP ksp, PetscInt restart)
276: {

280:   PetscTryMethod(ksp,"KSPGCRSetRestart_C",(KSP,PetscInt),(ksp,restart));
281:   return(0);
282: }

286: PetscErrorCode  KSPBuildSolution_GCR(KSP ksp, Vec v, Vec *V)
287: {
289:   KSP_GCR         *ctx;
290:   Vec             x;
291: 
293:   x = ksp->vec_sol;
294:   ctx = (KSP_GCR *)ksp->data;
295:   if (v) {
296:     VecCopy( x, v );
297:     if (V) *V = v;
298:   } else if (V) {
299:     *V = ksp->vec_sol;
300:   }
301:   return(0);
302: }

306: PetscErrorCode  KSPBuildResidual_GCR(KSP ksp, Vec t, Vec v, Vec *V)
307: {
309:   KSP_GCR         *ctx;
310: 
312:   ctx = (KSP_GCR *)ksp->data;
313:   if (v) {
314:     VecCopy( ctx->R, v );
315:     if (V) *V = v;
316:   } else if (V) {
317:     *V = ctx->R;
318:   }
319:   return(0);
320: }

322: /*MC
323:      KSPGCR - Implements the preconditioned Generalized Conjugate Residual method.  
324:  
325:  
326:    Options Database Keys:
327: +   -ksp_gcr_restart <restart> - the number of stored vectors to orthogonalize against
328:  
329:    Level: beginner
330:  
331:     Notes: The GCR Krylov method supports non-symmetric matrices and permits the use of a preconditioner 
332:            which may vary from one iteration to the next. Users can can define a method to vary the 
333:            preconditioner between iterates via KSPGCRSetModifyPC().
334:            Restarts are solves with x0 not equal to zero. When a restart occurs, the initial starting 
335:            solution is given by the current estimate for x which was obtained by the last restart 
336:            iterations of the GCR algorithm.
337:            When using GCR, the solution and residual vector can be directly accessed at any iterate,
338:            with zero computational cost, via a call to KSPBuildSolution() and KSPBuildResidual() respectively.
339:            This implementation of GCR will only apply the stopping condition test whenever ksp->its > ksp->chknorm, 
340:            where ksp->chknorm is specified via the command line argument -ksp_check_norm_iteration or via 
341:            the function KSPSetCheckNormIteration().
342:            The method implemented requires the storage of 2 x restart + 1 vectors.

344:     Contributed by Dave May
345:  
346:     References:
347:            S. C. Eisenstat, H. C. Elman, and H. C. Schultz. Variational iterative methods for 
348:            non-symmetric systems of linear equations. SIAM J. Numer. Anal., 20, 345-357, 1983
349:  
350:  
351: .seealso:  KSPCreate(), KSPSetType(), KSPType (for list of available types), KSP,
352:            KSPGCRSetRestart(), KSPGCRSetModifyPC(), KSPGMRES, KSPFGMRES
353:  
354: M*/
358: PetscErrorCode KSPCreate_GCR(KSP ksp)
359: {
361:   KSP_GCR        *ctx;

364:   PetscNewLog(ksp,KSP_GCR,&ctx);
365:   ctx->restart                   = 30;
366:   ctx->n_restarts                = 0;
367:   ksp->data                      = (void*)ctx;
368:   ksp->pc_side                   = PC_RIGHT;
369: 
370:   ksp->ops->setup                = KSPSetUp_GCR;
371:   ksp->ops->solve                = KSPSolve_GCR;
372:   ksp->ops->destroy              = KSPDestroy_GCR;
373:   ksp->ops->view                 = KSPView_GCR;
374:   ksp->ops->setfromoptions       = KSPSetFromOptions_GCR;
375:   ksp->ops->buildsolution        = KSPBuildSolution_GCR;
376:   ksp->ops->buildresidual        = KSPBuildResidual_GCR;
377: 
378:   PetscObjectComposeFunctionDynamic(  (PetscObject)ksp, "KSPGCRSetRestart_C",
379:                                       "KSPGCRSetRestart_GCR",KSPGCRSetRestart_GCR );
380:   PetscObjectComposeFunctionDynamic((PetscObject)ksp,"KSPGCRSetModifyPC_C",
381:                                            "KSPGCRSetModifyPC_GCR",KSPGCRSetModifyPC_GCR);
382:   return(0);
383: }