Actual source code: snesgs.c

  1: #include <../src/snes/impls/gs/gsimpl.h>

  3: /*@
  4:    SNESNGSSetTolerances - Sets various parameters used in convergence tests.

  6:    Logically Collective on SNES

  8:    Input Parameters:
  9: +  snes - the SNES context
 10: .  abstol - absolute convergence tolerance
 11: .  rtol - relative convergence tolerance
 12: .  stol -  convergence tolerance in terms of the norm of the change in the solution between steps,  || delta x || < stol*|| x ||
 13: -  maxit - maximum number of iterations

 15:    Options Database Keys:
 16: +    -snes_ngs_atol <abstol> - Sets abstol
 17: .    -snes_ngs_rtol <rtol> - Sets rtol
 18: .    -snes_ngs_stol <stol> - Sets stol
 19: -    -snes_max_it <maxit> - Sets maxit

 21:    Level: intermediate

 23: .seealso: SNESSetTrustRegionTolerance()
 24: @*/
 25: PetscErrorCode  SNESNGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)
 26: {
 27:   SNES_NGS *gs = (SNES_NGS*)snes->data;


 32:   if (abstol != PETSC_DEFAULT) {
 33:     if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
 34:     gs->abstol = abstol;
 35:   }
 36:   if (rtol != PETSC_DEFAULT) {
 37:     if (rtol < 0.0 || 1.0 <= rtol) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Relative tolerance %g must be non-negative and less than 1.0",(double)rtol);
 38:     gs->rtol = rtol;
 39:   }
 40:   if (stol != PETSC_DEFAULT) {
 41:     if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
 42:     gs->stol = stol;
 43:   }
 44:   if (maxit != PETSC_DEFAULT) {
 45:     if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
 46:     gs->max_its = maxit;
 47:   }
 48:   return(0);
 49: }

 51: /*@
 52:    SNESNGSGetTolerances - Gets various parameters used in convergence tests.

 54:    Not Collective

 56:    Input Parameters:
 57: +  snes - the SNES context
 58: .  atol - absolute convergence tolerance
 59: .  rtol - relative convergence tolerance
 60: .  stol -  convergence tolerance in terms of the norm
 61:            of the change in the solution between steps
 62: -  maxit - maximum number of iterations

 64:    Notes:
 65:    The user can specify NULL for any parameter that is not needed.

 67:    Level: intermediate

 69: .seealso: SNESSetTolerances()
 70: @*/
 71: PetscErrorCode  SNESNGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit)
 72: {
 73:   SNES_NGS *gs = (SNES_NGS*)snes->data;

 77:   if (atol) *atol = gs->abstol;
 78:   if (rtol) *rtol = gs->rtol;
 79:   if (stol) *stol = gs->stol;
 80:   if (maxit) *maxit = gs->max_its;
 81:   return(0);
 82: }

 84: /*@
 85:    SNESNGSSetSweeps - Sets the number of sweeps of GS to use.

 87:    Input Parameters:
 88: +  snes   - the SNES context
 89: -  sweeps  - the number of sweeps of GS to perform.

 91:    Level: intermediate

 93: .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetNPC(), SNESNGSGetSweeps()
 94: @*/

 96: PetscErrorCode SNESNGSSetSweeps(SNES snes, PetscInt sweeps)
 97: {
 98:   SNES_NGS *gs = (SNES_NGS*)snes->data;

102:   gs->sweeps = sweeps;
103:   return(0);
104: }

106: /*@
107:    SNESNGSGetSweeps - Gets the number of sweeps GS will use.

109:    Input Parameters:
110: .  snes   - the SNES context

112:    Output Parameters:
113: .  sweeps  - the number of sweeps of GS to perform.

115:    Level: intermediate

117: .seealso: SNESSetNGS(), SNESGetNGS(), SNESSetNPC(), SNESNGSSetSweeps()
118: @*/
119: PetscErrorCode SNESNGSGetSweeps(SNES snes, PetscInt * sweeps)
120: {
121:   SNES_NGS *gs = (SNES_NGS*)snes->data;

125:   *sweeps = gs->sweeps;
126:   return(0);
127: }

129: PetscErrorCode SNESReset_NGS(SNES snes)
130: {
131:   SNES_NGS       *gs = (SNES_NGS*)snes->data;

135:   ISColoringDestroy(&gs->coloring);
136:   return(0);
137: }

139: PetscErrorCode SNESDestroy_NGS(SNES snes)
140: {

144:   SNESReset_NGS(snes);
145:   PetscFree(snes->data);
146:   return(0);
147: }

149: PetscErrorCode SNESSetUp_NGS(SNES snes)
150: {
152:   PetscErrorCode (*f)(SNES,Vec,Vec,void*);

155:   SNESGetNGS(snes,&f,NULL);
156:   if (!f) {
157:     SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);
158:   }
159:   return(0);
160: }

162: PetscErrorCode SNESSetFromOptions_NGS(PetscOptionItems *PetscOptionsObject,SNES snes)
163: {
164:   SNES_NGS       *gs = (SNES_NGS*)snes->data;
166:   PetscInt       sweeps,max_its=PETSC_DEFAULT;
167:   PetscReal      rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT;
168:   PetscBool      flg,flg1,flg2,flg3;

171:   PetscOptionsHead(PetscOptionsObject,"SNES GS options");
172:   /* GS Options */
173:   PetscOptionsInt("-snes_ngs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);
174:   if (flg) {
175:     SNESNGSSetSweeps(snes,sweeps);
176:   }
177:   PetscOptionsReal("-snes_ngs_atol","Absolute residual tolerance for GS iteration","SNESComputeGS",gs->abstol,&atol,&flg);
178:   PetscOptionsReal("-snes_ngs_rtol","Relative residual tolerance for GS iteration","SNESComputeGS",gs->rtol,&rtol,&flg1);
179:   PetscOptionsReal("-snes_ngs_stol","Absolute update tolerance for GS iteration","SNESComputeGS",gs->stol,&stol,&flg2);
180:   PetscOptionsInt("-snes_ngs_max_it","Maximum number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3);
181:   if (flg || flg1 || flg2 || flg3) {
182:     SNESNGSSetTolerances(snes,atol,rtol,stol,max_its);
183:   }
184:   flg  = PETSC_FALSE;
185:   PetscOptionsBool("-snes_ngs_secant","Use finite difference secant approximation with coloring","",flg,&flg,NULL);
186:   if (flg) {
187:     SNESSetNGS(snes,SNESComputeNGSDefaultSecant,NULL);
188:     PetscInfo(snes,"Setting default finite difference secant approximation with coloring\n");
189:   }
190:   PetscOptionsReal("-snes_ngs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL);
191:   PetscOptionsBool("-snes_ngs_secant_mat_coloring","Use the graph coloring of the Jacobian for the secant GS","",gs->secant_mat,&gs->secant_mat,&flg);

193:   PetscOptionsTail();
194:   return(0);
195: }

197: PetscErrorCode SNESView_NGS(SNES snes, PetscViewer viewer)
198: {
200:   PetscErrorCode (*f)(SNES,Vec,Vec,void*);
201:   SNES_NGS       *gs = (SNES_NGS*)snes->data;
202:   PetscBool      iascii;

205:   PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);
206:   if (iascii) {
207:     DMSNESGetNGS(snes->dm,&f,NULL);
208:     if (f == SNESComputeNGSDefaultSecant) {
209:       PetscViewerASCIIPrintf(viewer,"  Use finite difference secant approximation with coloring with h = %g \n",(double)gs->h);
210:     }
211:   }
212:   return(0);
213: }

215: PetscErrorCode SNESSolve_NGS(SNES snes)
216: {
217:   Vec              F;
218:   Vec              X;
219:   Vec              B;
220:   PetscInt         i;
221:   PetscReal        fnorm;
222:   PetscErrorCode   ierr;
223:   SNESNormSchedule normschedule;


227:   if (snes->xl || snes->xu || snes->ops->computevariablebounds) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);

229:   PetscCitationsRegister(SNESCitation,&SNEScite);
230:   X = snes->vec_sol;
231:   F = snes->vec_func;
232:   B = snes->vec_rhs;

234:   PetscObjectSAWsTakeAccess((PetscObject)snes);
235:   snes->iter   = 0;
236:   snes->norm   = 0.;
237:   PetscObjectSAWsGrantAccess((PetscObject)snes);
238:   snes->reason = SNES_CONVERGED_ITERATING;

240:   SNESGetNormSchedule(snes, &normschedule);
241:   if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
242:     /* compute the initial function and preconditioned update delX */
243:     if (!snes->vec_func_init_set) {
244:       SNESComputeFunction(snes,X,F);
245:     } else snes->vec_func_init_set = PETSC_FALSE;

247:     VecNorm(F, NORM_2, &fnorm); /* fnorm <- ||F||  */
248:     SNESCheckFunctionNorm(snes,fnorm);
249:     PetscObjectSAWsTakeAccess((PetscObject)snes);
250:     snes->iter = 0;
251:     snes->norm = fnorm;
252:     PetscObjectSAWsGrantAccess((PetscObject)snes);
253:     SNESLogConvergenceHistory(snes,snes->norm,0);
254:     SNESMonitor(snes,0,snes->norm);

256:     /* test convergence */
257:     (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
258:     if (snes->reason) return(0);
259:   } else {
260:     PetscObjectSAWsGrantAccess((PetscObject)snes);
261:     SNESLogConvergenceHistory(snes,snes->norm,0);
262:   }

264:   /* Call general purpose update function */
265:   if (snes->ops->update) {
266:     (*snes->ops->update)(snes, snes->iter);
267:   }

269:   for (i = 0; i < snes->max_its; i++) {
270:     SNESComputeNGS(snes, B, X);
271:     /* only compute norms if requested or about to exit due to maximum iterations */
272:     if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
273:       SNESComputeFunction(snes,X,F);
274:       VecNorm(F, NORM_2, &fnorm); /* fnorm <- ||F||  */
275:       SNESCheckFunctionNorm(snes,fnorm);
276:       /* Monitor convergence */
277:       PetscObjectSAWsTakeAccess((PetscObject)snes);
278:       snes->iter = i+1;
279:       snes->norm = fnorm;
280:       PetscObjectSAWsGrantAccess((PetscObject)snes);
281:       SNESLogConvergenceHistory(snes,snes->norm,0);
282:       SNESMonitor(snes,snes->iter,snes->norm);
283:     }
284:     /* Test for convergence */
285:     if (normschedule == SNES_NORM_ALWAYS) {(*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);}
286:     if (snes->reason) return(0);
287:     /* Call general purpose update function */
288:     if (snes->ops->update) {
289:       (*snes->ops->update)(snes, snes->iter);
290:     }
291:   }
292:   if (normschedule == SNES_NORM_ALWAYS) {
293:     if (i == snes->max_its) {
294:       PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);
295:       if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
296:     }
297:   } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
298:   return(0);
299: }

301: /*MC
302:   SNESNGS - Either calls the user-provided solution routine provided with SNESSetNGS() or does a finite difference secant approximation
303:             using coloring.

305:    Level: advanced

307:   Options Database:
308: +   -snes_ngs_sweeps <n> - Number of sweeps of GS to apply
309: .    -snes_ngs_atol <atol> - Absolute residual tolerance for GS iteration
310: .    -snes_ngs_rtol <rtol> - Relative residual tolerance for GS iteration
311: .    -snes_ngs_stol <stol> - Absolute update tolerance for GS iteration
312: .    -snes_ngs_max_it <maxit> - Maximum number of sweeps of GS to apply
313: .    -snes_ngs_secant - Use pointwise secant local Jacobian approximation with coloring instead of user provided Gauss-Seidel routine, this is
314:                         used by default if no user provided Gauss-Seidel routine is available. Requires either that a DM that can compute a coloring
315:                         is available or a Jacobian sparse matrix is provided (from which to get the coloring).
316: .    -snes_ngs_secant_h <h> - Differencing parameter for secant approximation
317: .    -snes_ngs_secant_mat_coloring - Use the graph coloring of the Jacobian for the secant GS even if a DM is available.
318: -    -snes_norm_schedule <none, always, initialonly, finalonly, initialfinalonly> - how often the residual norms are computed

320:   Notes:
321:   the Gauss-Seidel smoother is inherited through composition.  If a solver has been created with SNESGetNPC(), it will have
322:   its parent's Gauss-Seidel routine associated with it.

324:   By default this routine computes the solution norm at each iteration, this can be time consuming, you can turn this off with SNESSetNormSchedule() or -snes_norm_schedule
325:    References:
326: .  1. - Peter R. Brune, Matthew G. Knepley, Barry F. Smith, and Xuemin Tu, "Composing Scalable Nonlinear Algebraic Solvers",
327:    SIAM Review, 57(4), 2015

329: .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetNGS(), SNESType (for list of available types), SNESNGSSetSweeps(), SNESNGSSetTolerances(),
330:           SNESSetNormSchedule()
331: M*/

333: PETSC_EXTERN PetscErrorCode SNESCreate_NGS(SNES snes)
334: {
335:   SNES_NGS        *gs;

339:   snes->ops->destroy        = SNESDestroy_NGS;
340:   snes->ops->setup          = SNESSetUp_NGS;
341:   snes->ops->setfromoptions = SNESSetFromOptions_NGS;
342:   snes->ops->view           = SNESView_NGS;
343:   snes->ops->solve          = SNESSolve_NGS;
344:   snes->ops->reset          = SNESReset_NGS;

346:   snes->usesksp = PETSC_FALSE;
347:   snes->usesnpc = PETSC_FALSE;

349:   snes->alwayscomputesfinalresidual = PETSC_FALSE;

351:   if (!snes->tolerancesset) {
352:     snes->max_its   = 10000;
353:     snes->max_funcs = 10000;
354:   }

356:   PetscNewLog(snes,&gs);

358:   gs->sweeps  = 1;
359:   gs->rtol    = 1e-5;
360:   gs->abstol  = PETSC_MACHINE_EPSILON;
361:   gs->stol    = 1000*PETSC_MACHINE_EPSILON;
362:   gs->max_its = 50;
363:   gs->h       = PETSC_SQRT_MACHINE_EPSILON;

365:   snes->data = (void*) gs;
366:   return(0);
367: }