xref: /petsc/src/snes/impls/ms/ms.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 #include <petsc/private/snesimpl.h> /*I "petscsnes.h" I*/
2 
3 static SNESMSType SNESMSDefault = SNESMSM62;
4 static PetscBool  SNESMSRegisterAllCalled;
5 static PetscBool  SNESMSPackageInitialized;
6 
7 typedef struct _SNESMSTableau *SNESMSTableau;
8 struct _SNESMSTableau {
9   char      *name;
10   PetscInt   nstages;    /* Number of stages */
11   PetscInt   nregisters; /* Number of registers */
12   PetscReal  stability;  /* Scaled stability region */
13   PetscReal *gamma;      /* Coefficients of 3S* method */
14   PetscReal *delta;      /* Coefficients of 3S* method */
15   PetscReal *betasub;    /* Subdiagonal of beta in Shu-Osher form */
16 };
17 
18 typedef struct _SNESMSTableauLink *SNESMSTableauLink;
19 struct _SNESMSTableauLink {
20   struct _SNESMSTableau tab;
21   SNESMSTableauLink     next;
22 };
23 static SNESMSTableauLink SNESMSTableauList;
24 
25 typedef struct {
26   SNESMSTableau tableau; /* Tableau in low-storage form */
27   PetscReal     damping; /* Damping parameter, like length of (pseudo) time step */
28   PetscBool     norms;   /* Compute norms, usually only for monitoring purposes */
29 } SNES_MS;
30 
31 /*@C
32   SNESMSRegisterAll - Registers all of the multi-stage methods in SNESMS
33 
34   Not Collective, but should be called by all processes which will need the schemes to be registered
35 
36   Level: advanced
37 
38 .seealso: `SNESMSRegisterDestroy()`
39 @*/
40 PetscErrorCode SNESMSRegisterAll(void) {
41   PetscFunctionBegin;
42   if (SNESMSRegisterAllCalled) PetscFunctionReturn(0);
43   SNESMSRegisterAllCalled = PETSC_TRUE;
44 
45   {
46     PetscReal alpha[1] = {1.0};
47     PetscCall(SNESMSRegister(SNESMSEULER, 1, 1, 1.0, NULL, NULL, alpha));
48   }
49 
50   {
51     PetscReal gamma[3][6] = {
52       {0.0000000000000000E+00, -7.0304722367110606E-01, -1.9836719667506464E-01, -1.6023843981863788E+00, 9.4483822882855284E-02,  -1.4204296130641869E-01},
53       {1.0000000000000000E+00, 1.1111025767083920E+00,  5.6150921583923230E-01,  7.4151723494934041E-01,  3.1714538168600587E-01,  4.6479276238548706E-01 },
54       {0.0000000000000000E+00, 0.0000000000000000E+00,  0.0000000000000000E+00,  6.7968174970583317E-01,  -4.1755042846051737E-03, -1.9115668129923846E-01}
55     };
56     PetscReal delta[6]   = {1.0000000000000000E+00, 5.3275427433201750E-01, 6.0143544663985238E-01, 4.5874077053842177E-01, 2.7544386906104651E-01, 0.0000000000000000E+00};
57     PetscReal betasub[6] = {8.4753115429481929E-01, 7.4018896368655618E-01, 6.5963574086583309E-03, 4.6747795645517759E-01, 1.3314545813643919E-01, 5.3260800028018784E-01};
58     PetscCall(SNESMSRegister(SNESMSM62, 6, 3, 1.0, &gamma[0][0], delta, betasub));
59   }
60 
61   { /* Jameson (1983) */
62     PetscReal alpha[4] = {0.25, 0.5, 0.55, 1.0};
63     PetscCall(SNESMSRegister(SNESMSJAMESON83, 4, 1, 1.0, NULL, NULL, alpha));
64   }
65 
66   { /* Van Leer, Tai, and Powell (1989) 1 stage, order 1 */
67     PetscReal alpha[1] = {1.0};
68     PetscCall(SNESMSRegister(SNESMSVLTP11, 1, 1, 0.5, NULL, NULL, alpha));
69   }
70   { /* Van Leer, Tai, and Powell (1989) 2 stage, order 1 */
71     PetscReal alpha[2] = {0.3333, 1.0};
72     PetscCall(SNESMSRegister(SNESMSVLTP21, 2, 1, 1.0, NULL, NULL, alpha));
73   }
74   { /* Van Leer, Tai, and Powell (1989) 3 stage, order 1 */
75     PetscReal alpha[3] = {0.1481, 0.4000, 1.0};
76     PetscCall(SNESMSRegister(SNESMSVLTP31, 3, 1, 1.5, NULL, NULL, alpha));
77   }
78   { /* Van Leer, Tai, and Powell (1989) 4 stage, order 1 */
79     PetscReal alpha[4] = {0.0833, 0.2069, 0.4265, 1.0};
80     PetscCall(SNESMSRegister(SNESMSVLTP41, 4, 1, 2.0, NULL, NULL, alpha));
81   }
82   { /* Van Leer, Tai, and Powell (1989) 5 stage, order 1 */
83     PetscReal alpha[5] = {0.0533, 0.1263, 0.2375, 0.4414, 1.0};
84     PetscCall(SNESMSRegister(SNESMSVLTP51, 5, 1, 2.5, NULL, NULL, alpha));
85   }
86   { /* Van Leer, Tai, and Powell (1989) 6 stage, order 1 */
87     PetscReal alpha[6] = {0.0370, 0.0851, 0.1521, 0.2562, 0.4512, 1.0};
88     PetscCall(SNESMSRegister(SNESMSVLTP61, 6, 1, 3.0, NULL, NULL, alpha));
89   }
90   PetscFunctionReturn(0);
91 }
92 
93 /*@C
94    SNESMSRegisterDestroy - Frees the list of schemes that were registered by SNESMSRegister().
95 
96    Not Collective
97 
98    Level: advanced
99 
100 .seealso: `SNESMSRegister()`, `SNESMSRegisterAll()`
101 @*/
102 PetscErrorCode SNESMSRegisterDestroy(void) {
103   SNESMSTableauLink link;
104 
105   PetscFunctionBegin;
106   while ((link = SNESMSTableauList)) {
107     SNESMSTableau t   = &link->tab;
108     SNESMSTableauList = link->next;
109 
110     PetscCall(PetscFree(t->name));
111     PetscCall(PetscFree(t->gamma));
112     PetscCall(PetscFree(t->delta));
113     PetscCall(PetscFree(t->betasub));
114     PetscCall(PetscFree(link));
115   }
116   SNESMSRegisterAllCalled = PETSC_FALSE;
117   PetscFunctionReturn(0);
118 }
119 
120 /*@C
121   SNESMSInitializePackage - This function initializes everything in the SNESMS package. It is called
122   from SNESInitializePackage().
123 
124   Level: developer
125 
126 .seealso: `PetscInitialize()`
127 @*/
128 PetscErrorCode SNESMSInitializePackage(void) {
129   PetscFunctionBegin;
130   if (SNESMSPackageInitialized) PetscFunctionReturn(0);
131   SNESMSPackageInitialized = PETSC_TRUE;
132 
133   PetscCall(SNESMSRegisterAll());
134   PetscCall(PetscRegisterFinalize(SNESMSFinalizePackage));
135   PetscFunctionReturn(0);
136 }
137 
138 /*@C
139   SNESMSFinalizePackage - This function destroys everything in the SNESMS package. It is
140   called from PetscFinalize().
141 
142   Level: developer
143 
144 .seealso: `PetscFinalize()`
145 @*/
146 PetscErrorCode SNESMSFinalizePackage(void) {
147   PetscFunctionBegin;
148   SNESMSPackageInitialized = PETSC_FALSE;
149 
150   PetscCall(SNESMSRegisterDestroy());
151   PetscFunctionReturn(0);
152 }
153 
154 /*@C
155    SNESMSRegister - register a multistage scheme
156 
157    Not Collective, but the same schemes should be registered on all processes on which they will be used
158 
159    Input Parameters:
160 +  name - identifier for method
161 .  nstages - number of stages
162 .  nregisters - number of registers used by low-storage implementation
163 .  stability - scaled stability region
164 .  gamma - coefficients, see Ketcheson's paper
165 .  delta - coefficients, see Ketcheson's paper
166 -  betasub - subdiagonal of Shu-Osher form
167 
168    Notes:
169    The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.
170 
171    Many multistage schemes are of the form
172    $ X_0 = X^{(n)}
173    $ X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s
174    $ X^{(n+1)} = X_s
175    These methods can be registered with
176 .vb
177    SNESMSRegister("name",s,1,stability,NULL,NULL,alpha);
178 .ve
179 
180    Level: advanced
181 
182 .seealso: `SNESMS`
183 @*/
184 PetscErrorCode SNESMSRegister(SNESMSType name, PetscInt nstages, PetscInt nregisters, PetscReal stability, const PetscReal gamma[], const PetscReal delta[], const PetscReal betasub[]) {
185   SNESMSTableauLink link;
186   SNESMSTableau     t;
187 
188   PetscFunctionBegin;
189   PetscValidCharPointer(name, 1);
190   PetscCheck(nstages >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have at least one stage");
191   if (gamma || delta) {
192     PetscCheck(nregisters == 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 3-register form");
193     PetscValidRealPointer(gamma, 5);
194     PetscValidRealPointer(delta, 6);
195   } else {
196     PetscCheck(nregisters == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 1-register form");
197   }
198   PetscValidRealPointer(betasub, 7);
199 
200   PetscCall(SNESMSInitializePackage());
201   PetscCall(PetscNew(&link));
202   t = &link->tab;
203   PetscCall(PetscStrallocpy(name, &t->name));
204   t->nstages    = nstages;
205   t->nregisters = nregisters;
206   t->stability  = stability;
207 
208   if (gamma && delta) {
209     PetscCall(PetscMalloc1(nstages * nregisters, &t->gamma));
210     PetscCall(PetscMalloc1(nstages, &t->delta));
211     PetscCall(PetscArraycpy(t->gamma, gamma, nstages * nregisters));
212     PetscCall(PetscArraycpy(t->delta, delta, nstages));
213   }
214   PetscCall(PetscMalloc1(nstages, &t->betasub));
215   PetscCall(PetscArraycpy(t->betasub, betasub, nstages));
216 
217   link->next        = SNESMSTableauList;
218   SNESMSTableauList = link;
219   PetscFunctionReturn(0);
220 }
221 
222 /*
223   X - initial state, updated in-place.
224   F - residual, computed at the initial X on input
225 */
226 static PetscErrorCode SNESMSStep_3Sstar(SNES snes, Vec X, Vec F) {
227   SNES_MS         *ms    = (SNES_MS *)snes->data;
228   SNESMSTableau    t     = ms->tableau;
229   const PetscReal *gamma = t->gamma, *delta = t->delta, *betasub = t->betasub;
230   Vec              S1, S2, S3, Y;
231   PetscInt         i, nstages = t->nstages;
232 
233   PetscFunctionBegin;
234   Y  = snes->work[0];
235   S1 = X;
236   S2 = snes->work[1];
237   S3 = snes->work[2];
238   PetscCall(VecZeroEntries(S2));
239   PetscCall(VecCopy(X, S3));
240   for (i = 0; i < nstages; i++) {
241     Vec         Ss[4];
242     PetscScalar scoeff[4];
243 
244     Ss[0] = S1;
245     Ss[1] = S2;
246     Ss[2] = S3;
247     Ss[3] = Y;
248 
249     scoeff[0] = gamma[0 * nstages + i] - 1;
250     scoeff[1] = gamma[1 * nstages + i];
251     scoeff[2] = gamma[2 * nstages + i];
252     scoeff[3] = -betasub[i] * ms->damping;
253 
254     PetscCall(VecAXPY(S2, delta[i], S1));
255     if (i > 0) { PetscCall(SNESComputeFunction(snes, S1, F)); }
256     PetscCall(KSPSolve(snes->ksp, F, Y));
257     PetscCall(VecMAXPY(S1, 4, scoeff, Ss));
258   }
259   PetscFunctionReturn(0);
260 }
261 
262 /*
263   X - initial state, updated in-place.
264   F - residual, computed at the initial X on input
265 */
266 static PetscErrorCode SNESMSStep_Basic(SNES snes, Vec X, Vec F) {
267   SNES_MS         *ms    = (SNES_MS *)snes->data;
268   SNESMSTableau    tab   = ms->tableau;
269   const PetscReal *alpha = tab->betasub, h = ms->damping;
270   PetscInt         i, nstages              = tab->nstages;
271   Vec              X0 = snes->work[0];
272 
273   PetscFunctionBegin;
274   PetscCall(VecCopy(X, X0));
275   for (i = 0; i < nstages; i++) {
276     if (i > 0) { PetscCall(SNESComputeFunction(snes, X, F)); }
277     PetscCall(KSPSolve(snes->ksp, F, X));
278     PetscCall(VecAYPX(X, -alpha[i] * h, X0));
279   }
280   PetscFunctionReturn(0);
281 }
282 
283 static PetscErrorCode SNESMSStep_Step(SNES snes, Vec X, Vec F) {
284   SNES_MS      *ms  = (SNES_MS *)snes->data;
285   SNESMSTableau tab = ms->tableau;
286 
287   PetscFunctionBegin;
288   if (tab->gamma && tab->delta) {
289     PetscCall(SNESMSStep_3Sstar(snes, X, F));
290   } else {
291     PetscCall(SNESMSStep_Basic(snes, X, F));
292   }
293   PetscFunctionReturn(0);
294 }
295 
296 static PetscErrorCode SNESMSStep_Norms(SNES snes, PetscInt iter, Vec F) {
297   SNES_MS  *ms = (SNES_MS *)snes->data;
298   PetscReal fnorm;
299 
300   PetscFunctionBegin;
301   if (ms->norms) {
302     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
303     SNESCheckFunctionNorm(snes, fnorm);
304     /* Monitor convergence */
305     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
306     snes->iter = iter;
307     snes->norm = fnorm;
308     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
309     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
310     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
311     /* Test for convergence */
312     PetscUseTypeMethod(snes, converged, snes->iter, 0.0, 0.0, fnorm, &snes->reason, snes->cnvP);
313   } else if (iter > 0) {
314     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
315     snes->iter = iter;
316     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
317   }
318   PetscFunctionReturn(0);
319 }
320 
321 static PetscErrorCode SNESSolve_MS(SNES snes) {
322   SNES_MS *ms = (SNES_MS *)snes->data;
323   Vec      X = snes->vec_sol, F = snes->vec_func;
324   PetscInt i;
325 
326   PetscFunctionBegin;
327   PetscCheck(!snes->xl && !snes->xu && !snes->ops->computevariablebounds, PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_WRONGSTATE, "SNES solver %s does not support bounds", ((PetscObject)snes)->type_name);
328   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));
329 
330   snes->reason = SNES_CONVERGED_ITERATING;
331   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
332   snes->iter = 0;
333   snes->norm = 0;
334   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
335 
336   if (!snes->vec_func_init_set) {
337     PetscCall(SNESComputeFunction(snes, X, F));
338   } else snes->vec_func_init_set = PETSC_FALSE;
339 
340   PetscCall(SNESMSStep_Norms(snes, 0, F));
341   if (snes->reason) PetscFunctionReturn(0);
342 
343   for (i = 0; i < snes->max_its; i++) {
344     /* Call general purpose update function */
345     PetscTryTypeMethod(snes, update, snes->iter);
346 
347     if (i == 0 && snes->jacobian) {
348       /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */
349       PetscCall(SNESComputeJacobian(snes, snes->vec_sol, snes->jacobian, snes->jacobian_pre));
350       SNESCheckJacobianDomainerror(snes);
351       PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
352     }
353 
354     PetscCall(SNESMSStep_Step(snes, X, F));
355 
356     if (i < snes->max_its - 1 || ms->norms) { PetscCall(SNESComputeFunction(snes, X, F)); }
357 
358     PetscCall(SNESMSStep_Norms(snes, i + 1, F));
359     if (snes->reason) PetscFunctionReturn(0);
360   }
361   if (!snes->reason) snes->reason = SNES_CONVERGED_ITS;
362   PetscFunctionReturn(0);
363 }
364 
365 static PetscErrorCode SNESSetUp_MS(SNES snes) {
366   SNES_MS      *ms    = (SNES_MS *)snes->data;
367   SNESMSTableau tab   = ms->tableau;
368   PetscInt      nwork = tab->nregisters;
369 
370   PetscFunctionBegin;
371   PetscCall(SNESSetWorkVecs(snes, nwork));
372   PetscCall(SNESSetUpMatrices(snes));
373   PetscFunctionReturn(0);
374 }
375 
376 static PetscErrorCode SNESReset_MS(SNES snes) {
377   PetscFunctionBegin;
378   PetscFunctionReturn(0);
379 }
380 
381 static PetscErrorCode SNESDestroy_MS(SNES snes) {
382   PetscFunctionBegin;
383   PetscCall(SNESReset_MS(snes));
384   PetscCall(PetscFree(snes->data));
385   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", NULL));
386   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", NULL));
387   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", NULL));
388   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", NULL));
389   PetscFunctionReturn(0);
390 }
391 
392 static PetscErrorCode SNESView_MS(SNES snes, PetscViewer viewer) {
393   PetscBool     iascii;
394   SNES_MS      *ms  = (SNES_MS *)snes->data;
395   SNESMSTableau tab = ms->tableau;
396 
397   PetscFunctionBegin;
398   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
399   if (iascii) { PetscCall(PetscViewerASCIIPrintf(viewer, "  multi-stage method type: %s\n", tab->name)); }
400   PetscFunctionReturn(0);
401 }
402 
403 static PetscErrorCode SNESSetFromOptions_MS(SNES snes, PetscOptionItems *PetscOptionsObject) {
404   SNES_MS *ms = (SNES_MS *)snes->data;
405 
406   PetscFunctionBegin;
407   PetscOptionsHeadBegin(PetscOptionsObject, "SNES MS options");
408   {
409     SNESMSTableauLink link;
410     PetscInt          count, choice;
411     PetscBool         flg;
412     const char      **namelist;
413     SNESMSType        mstype;
414     PetscReal         damping;
415 
416     PetscCall(SNESMSGetType(snes, &mstype));
417     for (link = SNESMSTableauList, count = 0; link; link = link->next, count++)
418       ;
419     PetscCall(PetscMalloc1(count, (char ***)&namelist));
420     for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
421     PetscCall(PetscOptionsEList("-snes_ms_type", "Multistage smoother type", "SNESMSSetType", (const char *const *)namelist, count, mstype, &choice, &flg));
422     if (flg) PetscCall(SNESMSSetType(snes, namelist[choice]));
423     PetscCall(PetscFree(namelist));
424     PetscCall(SNESMSGetDamping(snes, &damping));
425     PetscCall(PetscOptionsReal("-snes_ms_damping", "Damping for multistage method", "SNESMSSetDamping", damping, &damping, &flg));
426     if (flg) PetscCall(SNESMSSetDamping(snes, damping));
427     PetscCall(PetscOptionsBool("-snes_ms_norms", "Compute norms for monitoring", "none", ms->norms, &ms->norms, NULL));
428   }
429   PetscOptionsHeadEnd();
430   PetscFunctionReturn(0);
431 }
432 
433 static PetscErrorCode SNESMSGetType_MS(SNES snes, SNESMSType *mstype) {
434   SNES_MS      *ms  = (SNES_MS *)snes->data;
435   SNESMSTableau tab = ms->tableau;
436 
437   PetscFunctionBegin;
438   *mstype = tab->name;
439   PetscFunctionReturn(0);
440 }
441 
442 static PetscErrorCode SNESMSSetType_MS(SNES snes, SNESMSType mstype) {
443   SNES_MS          *ms = (SNES_MS *)snes->data;
444   SNESMSTableauLink link;
445   PetscBool         match;
446 
447   PetscFunctionBegin;
448   if (ms->tableau) {
449     PetscCall(PetscStrcmp(ms->tableau->name, mstype, &match));
450     if (match) PetscFunctionReturn(0);
451   }
452   for (link = SNESMSTableauList; link; link = link->next) {
453     PetscCall(PetscStrcmp(link->tab.name, mstype, &match));
454     if (match) {
455       if (snes->setupcalled) PetscCall(SNESReset_MS(snes));
456       ms->tableau = &link->tab;
457       if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes));
458       PetscFunctionReturn(0);
459     }
460   }
461   SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mstype);
462 }
463 
464 /*@C
465   SNESMSGetType - Get the type of multistage smoother
466 
467   Not collective
468 
469   Input Parameter:
470 .  snes - nonlinear solver context
471 
472   Output Parameter:
473 .  mstype - type of multistage method
474 
475   Level: beginner
476 
477 .seealso: `SNESMSSetType()`, `SNESMSType`, `SNESMS`
478 @*/
479 PetscErrorCode SNESMSGetType(SNES snes, SNESMSType *mstype) {
480   PetscFunctionBegin;
481   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
482   PetscValidPointer(mstype, 2);
483   PetscUseMethod(snes, "SNESMSGetType_C", (SNES, SNESMSType *), (snes, mstype));
484   PetscFunctionReturn(0);
485 }
486 
487 /*@C
488   SNESMSSetType - Set the type of multistage smoother
489 
490   Logically collective
491 
492   Input Parameters:
493 +  snes - nonlinear solver context
494 -  mstype - type of multistage method
495 
496   Level: beginner
497 
498 .seealso: `SNESMSGetType()`, `SNESMSType`, `SNESMS`
499 @*/
500 PetscErrorCode SNESMSSetType(SNES snes, SNESMSType mstype) {
501   PetscFunctionBegin;
502   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
503   PetscValidCharPointer(mstype, 2);
504   PetscTryMethod(snes, "SNESMSSetType_C", (SNES, SNESMSType), (snes, mstype));
505   PetscFunctionReturn(0);
506 }
507 
508 static PetscErrorCode SNESMSGetDamping_MS(SNES snes, PetscReal *damping) {
509   SNES_MS *ms = (SNES_MS *)snes->data;
510 
511   PetscFunctionBegin;
512   *damping = ms->damping;
513   PetscFunctionReturn(0);
514 }
515 
516 static PetscErrorCode SNESMSSetDamping_MS(SNES snes, PetscReal damping) {
517   SNES_MS *ms = (SNES_MS *)snes->data;
518 
519   PetscFunctionBegin;
520   ms->damping = damping;
521   PetscFunctionReturn(0);
522 }
523 
524 /*@C
525   SNESMSGetDamping - Get the damping parameter
526 
527   Not collective
528 
529   Input Parameter:
530 .  snes - nonlinear solver context
531 
532   Output Parameter:
533 .  damping - damping parameter
534 
535   Level: advanced
536 
537 .seealso: `SNESMSSetDamping()`, `SNESMS`
538 @*/
539 PetscErrorCode SNESMSGetDamping(SNES snes, PetscReal *damping) {
540   PetscFunctionBegin;
541   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
542   PetscValidRealPointer(damping, 2);
543   PetscUseMethod(snes, "SNESMSGetDamping_C", (SNES, PetscReal *), (snes, damping));
544   PetscFunctionReturn(0);
545 }
546 
547 /*@C
548   SNESMSSetDamping - Set the damping parameter
549 
550   Logically collective
551 
552   Input Parameters:
553 +  snes - nonlinear solver context
554 -  damping - damping parameter
555 
556   Level: advanced
557 
558 .seealso: `SNESMSGetDamping()`, `SNESMS`
559 @*/
560 PetscErrorCode SNESMSSetDamping(SNES snes, PetscReal damping) {
561   PetscFunctionBegin;
562   PetscValidHeaderSpecific(snes, SNES_CLASSID, 1);
563   PetscValidLogicalCollectiveReal(snes, damping, 2);
564   PetscTryMethod(snes, "SNESMSSetDamping_C", (SNES, PetscReal), (snes, damping));
565   PetscFunctionReturn(0);
566 }
567 
568 /* -------------------------------------------------------------------------- */
569 /*MC
570       SNESMS - multi-stage smoothers
571 
572       Options Database:
573 
574 +     -snes_ms_type - type of multi-stage smoother
575 -     -snes_ms_damping - damping for multi-stage method
576 
577       Notes:
578       These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for
579       FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev).
580 
581       Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds.
582 
583       The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with SNESMSRegister().
584 
585       References:
586 +     * - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006).
587 .     * - Jameson (1983) Solution of the Euler equations for two dimensional transonic flow by a multigrid method (https://doi.org/10.1016/0096-3003(83)90019-X).
588 .     * - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772).
589 -     * - Van Leer, Tai, and Powell (1989) Design of optimally smoothing multi-stage schemes for the Euler equations (https://doi.org/10.2514/6.1989-1933).
590 
591       Level: beginner
592 
593 .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESMS`, `SNESFAS`, `KSPCHEBYSHEV`
594 
595 M*/
596 PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes) {
597   SNES_MS *ms;
598 
599   PetscFunctionBegin;
600   PetscCall(SNESMSInitializePackage());
601 
602   snes->ops->setup          = SNESSetUp_MS;
603   snes->ops->solve          = SNESSolve_MS;
604   snes->ops->destroy        = SNESDestroy_MS;
605   snes->ops->setfromoptions = SNESSetFromOptions_MS;
606   snes->ops->view           = SNESView_MS;
607   snes->ops->reset          = SNESReset_MS;
608 
609   snes->usesnpc = PETSC_FALSE;
610   snes->usesksp = PETSC_TRUE;
611 
612   snes->alwayscomputesfinalresidual = PETSC_FALSE;
613 
614   PetscCall(PetscNewLog(snes, &ms));
615   snes->data  = (void *)ms;
616   ms->damping = 0.9;
617   ms->norms   = PETSC_FALSE;
618 
619   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", SNESMSGetType_MS));
620   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", SNESMSSetType_MS));
621   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", SNESMSGetDamping_MS));
622   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", SNESMSSetDamping_MS));
623 
624   PetscCall(SNESMSSetType(snes, SNESMSDefault));
625   PetscFunctionReturn(0);
626 }
627