Actual source code: snesgs.c
petsc-dev 2014-02-02
1: #include <../src/snes/impls/gs/gsimpl.h> /*I "petscsnes.h" I*/
5: /*@
6: SNESGSSetTolerances - Sets various parameters used in convergence tests.
8: Logically Collective on SNES
10: Input Parameters:
11: + snes - the SNES context
12: . abstol - absolute convergence tolerance
13: . rtol - relative convergence tolerance
14: . stol - convergence tolerance in terms of the norm of the change in the solution between steps, || delta x || < stol*|| x ||
15: - maxit - maximum number of iterations
18: Options Database Keys:
19: + -snes_gs_atol <abstol> - Sets abstol
20: . -snes_gs_rtol <rtol> - Sets rtol
21: . -snes_gs_stol <stol> - Sets stol
22: - -snes_max_it <maxit> - Sets maxit
24: Level: intermediate
26: .keywords: SNES, nonlinear, gauss-seidel, set, convergence, tolerances
28: .seealso: SNESSetTrustRegionTolerance()
29: @*/
30: PetscErrorCode SNESGSSetTolerances(SNES snes,PetscReal abstol,PetscReal rtol,PetscReal stol,PetscInt maxit)
31: {
32: SNES_GS *gs = (SNES_GS*)snes->data;
37: if (abstol != PETSC_DEFAULT) {
38: if (abstol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Absolute tolerance %g must be non-negative",(double)abstol);
39: gs->abstol = abstol;
40: }
41: if (rtol != PETSC_DEFAULT) {
42: 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);
43: gs->rtol = rtol;
44: }
45: if (stol != PETSC_DEFAULT) {
46: if (stol < 0.0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Step tolerance %g must be non-negative",(double)stol);
47: gs->stol = stol;
48: }
49: if (maxit != PETSC_DEFAULT) {
50: if (maxit < 0) SETERRQ1(PetscObjectComm((PetscObject)snes),PETSC_ERR_ARG_OUTOFRANGE,"Maximum number of iterations %D must be non-negative",maxit);
51: gs->max_its = maxit;
52: }
53: return(0);
54: }
58: /*@
59: SNESGSGetTolerances - Gets various parameters used in convergence tests.
61: Not Collective
63: Input Parameters:
64: + snes - the SNES context
65: . atol - absolute convergence tolerance
66: . rtol - relative convergence tolerance
67: . stol - convergence tolerance in terms of the norm
68: of the change in the solution between steps
69: - maxit - maximum number of iterations
71: Notes:
72: The user can specify NULL for any parameter that is not needed.
74: Level: intermediate
76: .keywords: SNES, nonlinear, get, convergence, tolerances
78: .seealso: SNESSetTolerances()
79: @*/
80: PetscErrorCode SNESGSGetTolerances(SNES snes,PetscReal *atol,PetscReal *rtol,PetscReal *stol,PetscInt *maxit)
81: {
82: SNES_GS *gs = (SNES_GS*)snes->data;
86: if (atol) *atol = gs->abstol;
87: if (rtol) *rtol = gs->rtol;
88: if (stol) *stol = gs->stol;
89: if (maxit) *maxit = gs->max_its;
90: return(0);
91: }
95: /*@
96: SNESGSSetSweeps - Sets the number of sweeps of GS to use.
98: Input Parameters:
99: + snes - the SNES context
100: - sweeps - the number of sweeps of GS to perform.
102: Level: intermediate
104: .keywords: SNES, nonlinear, set, Gauss-Siedel
106: .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSGetSweeps()
107: @*/
109: PetscErrorCode SNESGSSetSweeps(SNES snes, PetscInt sweeps)
110: {
111: SNES_GS *gs = (SNES_GS*)snes->data;
115: gs->sweeps = sweeps;
116: return(0);
117: }
121: /*@
122: SNESGSGetSweeps - Gets the number of sweeps GS will use.
124: Input Parameters:
125: . snes - the SNES context
127: Output Parameters:
128: . sweeps - the number of sweeps of GS to perform.
130: Level: intermediate
132: .keywords: SNES, nonlinear, set, Gauss-Siedel
134: .seealso: SNESSetGS(), SNESGetGS(), SNESSetPC(), SNESGSSetSweeps()
135: @*/
136: PetscErrorCode SNESGSGetSweeps(SNES snes, PetscInt * sweeps)
137: {
138: SNES_GS *gs = (SNES_GS*)snes->data;
142: *sweeps = gs->sweeps;
143: return(0);
144: }
149: PetscErrorCode SNESDefaultApplyGS(SNES snes,Vec X,Vec F,void *ctx)
150: {
152: /* see if there's a coloring on the DM */
153: return(0);
154: }
158: PetscErrorCode SNESReset_GS(SNES snes)
159: {
161: return(0);
162: }
166: PetscErrorCode SNESDestroy_GS(SNES snes)
167: {
171: SNESReset_GS(snes);
172: PetscFree(snes->data);
173: return(0);
174: }
178: PetscErrorCode SNESSetUp_GS(SNES snes)
179: {
181: PetscErrorCode (*f)(SNES,Vec,Vec,void*);
184: SNESGetGS(snes,&f,NULL);
185: if (!f) {
186: SNESSetGS(snes,SNESComputeGSDefaultSecant,NULL);
187: }
188: return(0);
189: }
193: PetscErrorCode SNESSetFromOptions_GS(SNES snes)
194: {
195: SNES_GS *gs = (SNES_GS*)snes->data;
197: PetscInt sweeps,max_its=PETSC_DEFAULT;
198: PetscReal rtol=PETSC_DEFAULT,atol=PETSC_DEFAULT,stol=PETSC_DEFAULT;
199: PetscBool flg,flg1,flg2,flg3;
202: PetscOptionsHead("SNES GS options");
203: /* GS Options */
204: PetscOptionsInt("-snes_gs_sweeps","Number of sweeps of GS to apply","SNESComputeGS",gs->sweeps,&sweeps,&flg);
205: if (flg) {
206: SNESGSSetSweeps(snes,sweeps);
207: }
208: PetscOptionsReal("-snes_gs_atol","Number of sweeps of GS to apply","SNESComputeGS",gs->abstol,&atol,&flg);
209: PetscOptionsReal("-snes_gs_rtol","Number of sweeps of GS to apply","SNESComputeGS",gs->rtol,&rtol,&flg1);
210: PetscOptionsReal("-snes_gs_stol","Number of sweeps of GS to apply","SNESComputeGS",gs->stol,&stol,&flg2);
211: PetscOptionsInt("-snes_gs_max_it","Number of sweeps of GS to apply","SNESComputeGS",gs->max_its,&max_its,&flg3);
212: if (flg || flg1 || flg2 || flg3) {
213: SNESGSSetTolerances(snes,atol,rtol,stol,max_its);
214: }
215: flg = PETSC_FALSE;
216: PetscOptionsBool("-snes_gs_secant","Use pointwise secant local Jacobian approximation","",flg,&flg,NULL);
217: if (flg) {
218: SNESSetGS(snes,SNESComputeGSDefaultSecant,NULL);
219: PetscInfo(snes,"Setting default finite difference coloring Jacobian matrix\n");
220: }
221: PetscOptionsReal("-snes_gs_secant_h","Differencing parameter for secant search","",gs->h,&gs->h,NULL);
222: PetscOptionsBool("-snes_gs_secant_mat_coloring","Use the Jacobian coloring for the secant GS","",gs->secant_mat,&gs->secant_mat,&flg);
224: PetscOptionsTail();
225: return(0);
226: }
230: PetscErrorCode SNESView_GS(SNES snes, PetscViewer viewer)
231: {
233: return(0);
234: }
238: PetscErrorCode SNESSolve_GS(SNES snes)
239: {
240: Vec F;
241: Vec X;
242: Vec B;
243: PetscInt i;
244: PetscReal fnorm;
245: PetscErrorCode ierr;
246: SNESNormSchedule normschedule;
249: X = snes->vec_sol;
250: F = snes->vec_func;
251: B = snes->vec_rhs;
253: PetscObjectSAWsTakeAccess((PetscObject)snes);
254: snes->iter = 0;
255: snes->norm = 0.;
256: PetscObjectSAWsGrantAccess((PetscObject)snes);
257: snes->reason = SNES_CONVERGED_ITERATING;
259: SNESGetNormSchedule(snes, &normschedule);
260: if (normschedule == SNES_NORM_ALWAYS || normschedule == SNES_NORM_INITIAL_ONLY || normschedule == SNES_NORM_INITIAL_FINAL_ONLY) {
261: /* compute the initial function and preconditioned update delX */
262: if (!snes->vec_func_init_set) {
263: SNESComputeFunction(snes,X,F);
264: if (snes->domainerror) {
265: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
266: return(0);
267: }
268: } else snes->vec_func_init_set = PETSC_FALSE;
270: VecNorm(F, NORM_2, &fnorm); /* fnorm <- ||F|| */
271: if (PetscIsInfOrNanReal(fnorm)) {
272: snes->reason = SNES_DIVERGED_FNORM_NAN;
273: return(0);
274: }
275: PetscObjectSAWsTakeAccess((PetscObject)snes);
276: snes->iter = 0;
277: snes->norm = fnorm;
278: PetscObjectSAWsGrantAccess((PetscObject)snes);
279: SNESLogConvergenceHistory(snes,snes->norm,0);
280: SNESMonitor(snes,0,snes->norm);
282: /* test convergence */
283: (*snes->ops->converged)(snes,0,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
284: if (snes->reason) return(0);
285: } else {
286: PetscObjectSAWsGrantAccess((PetscObject)snes);
287: SNESLogConvergenceHistory(snes,snes->norm,0);
288: SNESMonitor(snes,0,snes->norm);
289: }
291: /* Call general purpose update function */
292: if (snes->ops->update) {
293: (*snes->ops->update)(snes, snes->iter);
294: }
296: for (i = 0; i < snes->max_its; i++) {
297: SNESComputeGS(snes, B, X);
298: /* only compute norms if requested or about to exit due to maximum iterations */
299: if (normschedule == SNES_NORM_ALWAYS || ((i == snes->max_its - 1) && (normschedule == SNES_NORM_INITIAL_FINAL_ONLY || normschedule == SNES_NORM_FINAL_ONLY))) {
300: SNESComputeFunction(snes,X,F);
301: if (snes->domainerror) {
302: snes->reason = SNES_DIVERGED_FUNCTION_DOMAIN;
303: return(0);
304: }
305: VecNorm(F, NORM_2, &fnorm); /* fnorm <- ||F|| */
306: if (PetscIsInfOrNanReal(fnorm)) {
307: snes->reason = SNES_DIVERGED_FNORM_NAN;
308: return(0);
309: }
310: }
311: /* Monitor convergence */
312: PetscObjectSAWsTakeAccess((PetscObject)snes);
313: snes->iter = i+1;
314: snes->norm = fnorm;
315: PetscObjectSAWsGrantAccess((PetscObject)snes);
316: SNESLogConvergenceHistory(snes,snes->norm,0);
317: SNESMonitor(snes,snes->iter,snes->norm);
318: /* Test for convergence */
319: if (normschedule == SNES_NORM_ALWAYS) (*snes->ops->converged)(snes,snes->iter,0.0,0.0,fnorm,&snes->reason,snes->cnvP);
320: if (snes->reason) return(0);
321: /* Call general purpose update function */
322: if (snes->ops->update) {
323: (*snes->ops->update)(snes, snes->iter);
324: }
325: }
326: if (normschedule == SNES_NORM_ALWAYS) {
327: if (i == snes->max_its) {
328: PetscInfo1(snes,"Maximum number of iterations has been reached: %D\n",snes->max_its);
329: if (!snes->reason) snes->reason = SNES_DIVERGED_MAX_IT;
330: }
331: } else if (!snes->reason) snes->reason = SNES_CONVERGED_ITS; /* GS is meant to be used as a preconditioner */
332: return(0);
333: }
335: /*MC
336: SNESGS - Just calls the user-provided solution routine provided with SNESSetGS()
338: Level: advanced
340: Notes:
341: the Gauss-Seidel smoother is inherited through composition. If a solver has been created with SNESGetPC(), it will have
342: its parent's Gauss-Seidel routine associated with it.
344: .seealso: SNESCreate(), SNES, SNESSetType(), SNESSetGS(), SNESType (for list of available types)
345: M*/
349: PETSC_EXTERN PetscErrorCode SNESCreate_GS(SNES snes)
350: {
351: SNES_GS *gs;
355: snes->ops->destroy = SNESDestroy_GS;
356: snes->ops->setup = SNESSetUp_GS;
357: snes->ops->setfromoptions = SNESSetFromOptions_GS;
358: snes->ops->view = SNESView_GS;
359: snes->ops->solve = SNESSolve_GS;
360: snes->ops->reset = SNESReset_GS;
362: snes->usesksp = PETSC_FALSE;
363: snes->usespc = PETSC_FALSE;
365: if (!snes->tolerancesset) {
366: snes->max_its = 10000;
367: snes->max_funcs = 10000;
368: }
370: PetscNewLog(snes,&gs);
372: gs->sweeps = 1;
373: gs->rtol = 1e-5;
374: gs->abstol = 1e-15;
375: gs->stol = 1e-12;
376: gs->max_its = 50;
377: gs->h = 1e-8;
379: snes->data = (void*) gs;
380: return(0);
381: }