Actual source code: ms.c

  1: #include <petsc/private/snesimpl.h>

  3: static SNESMSType SNESMSDefault = SNESMSM62;
  4: static PetscBool  SNESMSRegisterAllCalled;
  5: static PetscBool  SNESMSPackageInitialized;

  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: };

 18: typedef struct _SNESMSTableauLink *SNESMSTableauLink;
 19: struct _SNESMSTableauLink {
 20:   struct _SNESMSTableau tab;
 21:   SNESMSTableauLink     next;
 22: };
 23: static SNESMSTableauLink SNESMSTableauList;

 25: typedef struct {
 26:   SNESMSTableau tableau; /* Tableau in low-storage form */
 27:   PetscReal     damping; /* Damping parameter, like length of (pseudo) time step */
 28: } SNES_MS;

 30: /*@C
 31:   SNESMSRegisterAll - Registers all of the multi-stage methods in `SNESMS`

 33:   Logically Collective

 35:   Level: developer

 37: .seealso: `SNES`, `SNESMS`, `SNESMSRegisterDestroy()`
 38: @*/
 39: PetscErrorCode SNESMSRegisterAll(void)
 40: {
 41:   PetscFunctionBegin;
 42:   if (SNESMSRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
 43:   SNESMSRegisterAllCalled = PETSC_TRUE;

 45:   {
 46:     PetscReal alpha[1] = {1.0};
 47:     PetscCall(SNESMSRegister(SNESMSEULER, 1, 1, 1.0, NULL, NULL, alpha));
 48:   }

 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:   }

 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:   }

 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(PETSC_SUCCESS);
 91: }

 93: /*@C
 94:   SNESMSRegisterDestroy - Frees the list of schemes that were registered by `SNESMSRegister()`.

 96:   Not Collective

 98:   Level: developer

100: .seealso: `SNESMSRegister()`, `SNESMSRegisterAll()`
101: @*/
102: PetscErrorCode SNESMSRegisterDestroy(void)
103: {
104:   SNESMSTableauLink link;

106:   PetscFunctionBegin;
107:   while ((link = SNESMSTableauList)) {
108:     SNESMSTableau t   = &link->tab;
109:     SNESMSTableauList = link->next;

111:     PetscCall(PetscFree(t->name));
112:     PetscCall(PetscFree(t->gamma));
113:     PetscCall(PetscFree(t->delta));
114:     PetscCall(PetscFree(t->betasub));
115:     PetscCall(PetscFree(link));
116:   }
117:   SNESMSRegisterAllCalled = PETSC_FALSE;
118:   PetscFunctionReturn(PETSC_SUCCESS);
119: }

121: /*@C
122:   SNESMSInitializePackage - This function initializes everything in the `SNESMS` package. It is called
123:   from `SNESInitializePackage()`.

125:   Level: developer

127: .seealso: `PetscInitialize()`
128: @*/
129: PetscErrorCode SNESMSInitializePackage(void)
130: {
131:   PetscFunctionBegin;
132:   if (SNESMSPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
133:   SNESMSPackageInitialized = PETSC_TRUE;

135:   PetscCall(SNESMSRegisterAll());
136:   PetscCall(PetscRegisterFinalize(SNESMSFinalizePackage));
137:   PetscFunctionReturn(PETSC_SUCCESS);
138: }

140: /*@C
141:   SNESMSFinalizePackage - This function destroys everything in the `SNESMS` package. It is
142:   called from `PetscFinalize()`.

144:   Level: developer

146: .seealso: `PetscFinalize()`
147: @*/
148: PetscErrorCode SNESMSFinalizePackage(void)
149: {
150:   PetscFunctionBegin;
151:   SNESMSPackageInitialized = PETSC_FALSE;

153:   PetscCall(SNESMSRegisterDestroy());
154:   PetscFunctionReturn(PETSC_SUCCESS);
155: }

157: /*@C
158:   SNESMSRegister - register a multistage scheme for `SNESMS`

160:   Logically Collective

162:   Input Parameters:
163: + name       - identifier for method
164: . nstages    - number of stages
165: . nregisters - number of registers used by low-storage implementation
166: . stability  - scaled stability region
167: . gamma      - coefficients, see Ketcheson's paper
168: . delta      - coefficients, see Ketcheson's paper
169: - betasub    - subdiagonal of Shu-Osher form

171:   Level: advanced

173:   Notes:
174:   The notation is described in Ketcheson (2010) Runge-Kutta methods with minimum storage implementations.

176:   Many multistage schemes are of the form
177:    $ X_0 = X^{(n)}
178:    $ X_k = X_0 + \alpha_k * F(X_{k-1}), k = 1,\ldots,s
179:    $ X^{(n+1)} = X_s
180:   These methods can be registered with
181: .vb
182:    SNESMSRegister("name",s,1,stability,NULL,NULL,alpha);
183: .ve

185: .seealso: `SNESMS`
186: @*/
187: PetscErrorCode SNESMSRegister(SNESMSType name, PetscInt nstages, PetscInt nregisters, PetscReal stability, const PetscReal gamma[], const PetscReal delta[], const PetscReal betasub[])
188: {
189:   SNESMSTableauLink link;
190:   SNESMSTableau     t;

192:   PetscFunctionBegin;
193:   PetscAssertPointer(name, 1);
194:   PetscCheck(nstages >= 1, PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "Must have at least one stage");
195:   if (gamma || delta) {
196:     PetscCheck(nregisters == 3, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 3-register form");
197:     PetscAssertPointer(gamma, 5);
198:     PetscAssertPointer(delta, 6);
199:   } else {
200:     PetscCheck(nregisters == 1, PETSC_COMM_SELF, PETSC_ERR_SUP, "Only support for methods written in 1-register form");
201:   }
202:   PetscAssertPointer(betasub, 7);

204:   PetscCall(SNESMSInitializePackage());
205:   PetscCall(PetscNew(&link));
206:   t = &link->tab;
207:   PetscCall(PetscStrallocpy(name, &t->name));
208:   t->nstages    = nstages;
209:   t->nregisters = nregisters;
210:   t->stability  = stability;

212:   if (gamma && delta) {
213:     PetscCall(PetscMalloc1(nstages * nregisters, &t->gamma));
214:     PetscCall(PetscMalloc1(nstages, &t->delta));
215:     PetscCall(PetscArraycpy(t->gamma, gamma, nstages * nregisters));
216:     PetscCall(PetscArraycpy(t->delta, delta, nstages));
217:   }
218:   PetscCall(PetscMalloc1(nstages, &t->betasub));
219:   PetscCall(PetscArraycpy(t->betasub, betasub, nstages));

221:   link->next        = SNESMSTableauList;
222:   SNESMSTableauList = link;
223:   PetscFunctionReturn(PETSC_SUCCESS);
224: }

226: /*
227:   X - initial state, updated in-place.
228:   F - residual, computed at the initial X on input
229: */
230: static PetscErrorCode SNESMSStep_3Sstar(SNES snes, Vec X, Vec F)
231: {
232:   SNES_MS         *ms      = (SNES_MS *)snes->data;
233:   SNESMSTableau    t       = ms->tableau;
234:   const PetscInt   nstages = t->nstages;
235:   const PetscReal *gamma = t->gamma, *delta = t->delta, *betasub = t->betasub;
236:   Vec              S1 = X, S2 = snes->work[1], S3 = snes->work[2], Y = snes->work[0], S1copy = snes->work[3];

238:   PetscFunctionBegin;
239:   PetscCall(VecZeroEntries(S2));
240:   PetscCall(VecCopy(X, S3));
241:   for (PetscInt i = 0; i < nstages; i++) {
242:     Vec               Ss[]     = {S1copy, S2, S3, Y};
243:     const PetscScalar scoeff[] = {gamma[0 * nstages + i] - 1, gamma[1 * nstages + i], gamma[2 * nstages + i], -betasub[i] * ms->damping};

245:     PetscCall(VecAXPY(S2, delta[i], S1));
246:     if (i > 0) PetscCall(SNESComputeFunction(snes, S1, F));
247:     PetscCall(KSPSolve(snes->ksp, F, Y));
248:     PetscCall(VecCopy(S1, S1copy));
249:     PetscCall(VecMAXPY(S1, 4, scoeff, Ss));
250:   }
251:   PetscFunctionReturn(PETSC_SUCCESS);
252: }

254: /*
255:   X - initial state, updated in-place.
256:   F - residual, computed at the initial X on input
257: */
258: static PetscErrorCode SNESMSStep_Basic(SNES snes, Vec X, Vec F)
259: {
260:   SNES_MS         *ms    = (SNES_MS *)snes->data;
261:   SNESMSTableau    tab   = ms->tableau;
262:   const PetscReal *alpha = tab->betasub, h = ms->damping;
263:   PetscInt         i, nstages              = tab->nstages;
264:   Vec              X0 = snes->work[0];

266:   PetscFunctionBegin;
267:   PetscCall(VecCopy(X, X0));
268:   for (i = 0; i < nstages; i++) {
269:     if (i > 0) PetscCall(SNESComputeFunction(snes, X, F));
270:     PetscCall(KSPSolve(snes->ksp, F, X));
271:     PetscCall(VecAYPX(X, -alpha[i] * h, X0));
272:   }
273:   PetscFunctionReturn(PETSC_SUCCESS);
274: }

276: static PetscErrorCode SNESMSStep_Step(SNES snes, Vec X, Vec F)
277: {
278:   SNES_MS      *ms  = (SNES_MS *)snes->data;
279:   SNESMSTableau tab = ms->tableau;

281:   PetscFunctionBegin;
282:   if (tab->gamma && tab->delta) {
283:     PetscCall(SNESMSStep_3Sstar(snes, X, F));
284:   } else {
285:     PetscCall(SNESMSStep_Basic(snes, X, F));
286:   }
287:   PetscFunctionReturn(PETSC_SUCCESS);
288: }

290: static PetscErrorCode SNESMSStep_Norms(SNES snes, PetscInt iter, Vec F)
291: {
292:   PetscReal fnorm;

294:   PetscFunctionBegin;
295:   if (SNESNeedNorm_Private(snes, iter)) {
296:     PetscCall(VecNorm(F, NORM_2, &fnorm)); /* fnorm <- ||F||  */
297:     SNESCheckFunctionNorm(snes, fnorm);
298:     /* Monitor convergence */
299:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
300:     snes->iter = iter;
301:     snes->norm = fnorm;
302:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
303:     PetscCall(SNESLogConvergenceHistory(snes, snes->norm, 0));
304:     /* Test for convergence */
305:     PetscCall(SNESConverged(snes, snes->iter, 0.0, 0.0, fnorm));
306:     PetscCall(SNESMonitor(snes, snes->iter, snes->norm));
307:   } else if (iter > 0) {
308:     PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
309:     snes->iter = iter;
310:     PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));
311:   }
312:   PetscFunctionReturn(PETSC_SUCCESS);
313: }

315: static PetscErrorCode SNESSolve_MS(SNES snes)
316: {
317:   Vec      X = snes->vec_sol, F = snes->vec_func;
318:   PetscInt i;

320:   PetscFunctionBegin;
321:   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);
322:   PetscCall(PetscCitationsRegister(SNESCitation, &SNEScite));

324:   snes->reason = SNES_CONVERGED_ITERATING;
325:   PetscCall(PetscObjectSAWsTakeAccess((PetscObject)snes));
326:   snes->iter = 0;
327:   snes->norm = 0;
328:   PetscCall(PetscObjectSAWsGrantAccess((PetscObject)snes));

330:   if (!snes->vec_func_init_set) {
331:     PetscCall(SNESComputeFunction(snes, X, F));
332:   } else snes->vec_func_init_set = PETSC_FALSE;

334:   PetscCall(SNESMSStep_Norms(snes, 0, F));
335:   if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);

337:   for (i = 0; i < snes->max_its; i++) {
338:     /* Call general purpose update function */
339:     PetscTryTypeMethod(snes, update, snes->iter);

341:     if (i == 0 && snes->jacobian) {
342:       /* This method does not require a Jacobian, but it is usually preconditioned by PBJacobi */
343:       PetscCall(SNESComputeJacobian(snes, snes->vec_sol, snes->jacobian, snes->jacobian_pre));
344:       SNESCheckJacobianDomainerror(snes);
345:       PetscCall(KSPSetOperators(snes->ksp, snes->jacobian, snes->jacobian_pre));
346:     }

348:     PetscCall(SNESMSStep_Step(snes, X, F));

350:     if (i < snes->max_its - 1 || SNESNeedNorm_Private(snes, i + 1)) PetscCall(SNESComputeFunction(snes, X, F));

352:     PetscCall(SNESMSStep_Norms(snes, i + 1, F));
353:     if (snes->reason) PetscFunctionReturn(PETSC_SUCCESS);
354:   }
355:   PetscFunctionReturn(PETSC_SUCCESS);
356: }

358: static PetscErrorCode SNESSetUp_MS(SNES snes)
359: {
360:   SNES_MS      *ms    = (SNES_MS *)snes->data;
361:   SNESMSTableau tab   = ms->tableau;
362:   PetscInt      nwork = tab->nregisters + 1; // +1 because VecMAXPY() in SNESMSStep_3Sstar()
363:                                              // needs an extra work vec

365:   PetscFunctionBegin;
366:   PetscCall(SNESSetWorkVecs(snes, nwork));
367:   PetscCall(SNESSetUpMatrices(snes));
368:   PetscFunctionReturn(PETSC_SUCCESS);
369: }

371: static PetscErrorCode SNESReset_MS(SNES snes)
372: {
373:   PetscFunctionBegin;
374:   PetscFunctionReturn(PETSC_SUCCESS);
375: }

377: static PetscErrorCode SNESDestroy_MS(SNES snes)
378: {
379:   PetscFunctionBegin;
380:   PetscCall(SNESReset_MS(snes));
381:   PetscCall(PetscFree(snes->data));
382:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", NULL));
383:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", NULL));
384:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", NULL));
385:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", NULL));
386:   PetscFunctionReturn(PETSC_SUCCESS);
387: }

389: static PetscErrorCode SNESView_MS(SNES snes, PetscViewer viewer)
390: {
391:   PetscBool     iascii;
392:   SNES_MS      *ms  = (SNES_MS *)snes->data;
393:   SNESMSTableau tab = ms->tableau;

395:   PetscFunctionBegin;
396:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
397:   if (iascii) PetscCall(PetscViewerASCIIPrintf(viewer, "  multi-stage method type: %s\n", tab->name));
398:   PetscFunctionReturn(PETSC_SUCCESS);
399: }

401: static PetscErrorCode SNESSetFromOptions_MS(SNES snes, PetscOptionItems *PetscOptionsObject)
402: {
403:   PetscFunctionBegin;
404:   PetscOptionsHeadBegin(PetscOptionsObject, "SNES MS options");
405:   {
406:     SNESMSTableauLink link;
407:     PetscInt          count, choice;
408:     PetscBool         flg;
409:     const char      **namelist;
410:     SNESMSType        mstype;
411:     PetscReal         damping;
412:     PetscBool         norms = PETSC_FALSE;

414:     PetscCall(SNESMSGetType(snes, &mstype));
415:     for (link = SNESMSTableauList, count = 0; link; link = link->next, count++)
416:       ;
417:     PetscCall(PetscMalloc1(count, (char ***)&namelist));
418:     for (link = SNESMSTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
419:     PetscCall(PetscOptionsEList("-snes_ms_type", "Multistage smoother type", "SNESMSSetType", (const char *const *)namelist, count, mstype, &choice, &flg));
420:     if (flg) PetscCall(SNESMSSetType(snes, namelist[choice]));
421:     PetscCall(PetscFree(namelist));
422:     PetscCall(SNESMSGetDamping(snes, &damping));
423:     PetscCall(PetscOptionsReal("-snes_ms_damping", "Damping for multistage method", "SNESMSSetDamping", damping, &damping, &flg));
424:     if (flg) PetscCall(SNESMSSetDamping(snes, damping));

426:     /* deprecated option */
427:     PetscCall(PetscOptionsDeprecated("-snes_ms_norms", NULL, "3.20", "Use -snes_norm_schedule always"));
428:     PetscCall(PetscOptionsBool("-snes_ms_norms", NULL, NULL, norms, &norms, NULL));
429:     if (norms) PetscCall(SNESSetNormSchedule(snes, SNES_NORM_ALWAYS));
430:   }
431:   PetscOptionsHeadEnd();
432:   PetscFunctionReturn(PETSC_SUCCESS);
433: }

435: static PetscErrorCode SNESMSGetType_MS(SNES snes, SNESMSType *mstype)
436: {
437:   SNES_MS      *ms  = (SNES_MS *)snes->data;
438:   SNESMSTableau tab = ms->tableau;

440:   PetscFunctionBegin;
441:   *mstype = tab->name;
442:   PetscFunctionReturn(PETSC_SUCCESS);
443: }

445: static PetscErrorCode SNESMSSetType_MS(SNES snes, SNESMSType mstype)
446: {
447:   SNES_MS          *ms = (SNES_MS *)snes->data;
448:   SNESMSTableauLink link;
449:   PetscBool         match;

451:   PetscFunctionBegin;
452:   if (ms->tableau) {
453:     PetscCall(PetscStrcmp(ms->tableau->name, mstype, &match));
454:     if (match) PetscFunctionReturn(PETSC_SUCCESS);
455:   }
456:   for (link = SNESMSTableauList; link; link = link->next) {
457:     PetscCall(PetscStrcmp(link->tab.name, mstype, &match));
458:     if (match) {
459:       if (snes->setupcalled) PetscCall(SNESReset_MS(snes));
460:       ms->tableau = &link->tab;
461:       if (snes->setupcalled) PetscCall(SNESSetUp_MS(snes));
462:       PetscFunctionReturn(PETSC_SUCCESS);
463:     }
464:   }
465:   SETERRQ(PetscObjectComm((PetscObject)snes), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", mstype);
466: }

468: /*@C
469:   SNESMSGetType - Get the type of multistage smoother `SNESMS`

471:   Not Collective

473:   Input Parameter:
474: . snes - nonlinear solver context

476:   Output Parameter:
477: . mstype - type of multistage method

479:   Level: advanced

481: .seealso: `SNESMS`, `SNESMSSetType()`, `SNESMSType`
482: @*/
483: PetscErrorCode SNESMSGetType(SNES snes, SNESMSType *mstype)
484: {
485:   PetscFunctionBegin;
487:   PetscAssertPointer(mstype, 2);
488:   PetscUseMethod(snes, "SNESMSGetType_C", (SNES, SNESMSType *), (snes, mstype));
489:   PetscFunctionReturn(PETSC_SUCCESS);
490: }

492: /*@C
493:   SNESMSSetType - Set the type of multistage smoother `SNESMS`

495:   Logically Collective

497:   Input Parameters:
498: + snes   - nonlinear solver context
499: - mstype - type of multistage method

501:   Level: advanced

503: .seealso: `SNESMS`, `SNESMSGetType()`, `SNESMSType`
504: @*/
505: PetscErrorCode SNESMSSetType(SNES snes, SNESMSType mstype)
506: {
507:   PetscFunctionBegin;
509:   PetscAssertPointer(mstype, 2);
510:   PetscTryMethod(snes, "SNESMSSetType_C", (SNES, SNESMSType), (snes, mstype));
511:   PetscFunctionReturn(PETSC_SUCCESS);
512: }

514: static PetscErrorCode SNESMSGetDamping_MS(SNES snes, PetscReal *damping)
515: {
516:   SNES_MS *ms = (SNES_MS *)snes->data;

518:   PetscFunctionBegin;
519:   *damping = ms->damping;
520:   PetscFunctionReturn(PETSC_SUCCESS);
521: }

523: static PetscErrorCode SNESMSSetDamping_MS(SNES snes, PetscReal damping)
524: {
525:   SNES_MS *ms = (SNES_MS *)snes->data;

527:   PetscFunctionBegin;
528:   ms->damping = damping;
529:   PetscFunctionReturn(PETSC_SUCCESS);
530: }

532: /*@C
533:   SNESMSGetDamping - Get the damping parameter of `SNESMS` multistage scheme

535:   Not Collective

537:   Input Parameter:
538: . snes - nonlinear solver context

540:   Output Parameter:
541: . damping - damping parameter

543:   Level: advanced

545: .seealso: `SNESMSSetDamping()`, `SNESMS`
546: @*/
547: PetscErrorCode SNESMSGetDamping(SNES snes, PetscReal *damping)
548: {
549:   PetscFunctionBegin;
551:   PetscAssertPointer(damping, 2);
552:   PetscUseMethod(snes, "SNESMSGetDamping_C", (SNES, PetscReal *), (snes, damping));
553:   PetscFunctionReturn(PETSC_SUCCESS);
554: }

556: /*@C
557:   SNESMSSetDamping - Set the damping parameter for a `SNESMS` multistage scheme

559:   Logically Collective

561:   Input Parameters:
562: + snes    - nonlinear solver context
563: - damping - damping parameter

565:   Level: advanced

567: .seealso: `SNESMSGetDamping()`, `SNESMS`
568: @*/
569: PetscErrorCode SNESMSSetDamping(SNES snes, PetscReal damping)
570: {
571:   PetscFunctionBegin;
574:   PetscTryMethod(snes, "SNESMSSetDamping_C", (SNES, PetscReal), (snes, damping));
575:   PetscFunctionReturn(PETSC_SUCCESS);
576: }

578: /*MC
579:       SNESMS - multi-stage smoothers

581:       Options Database Keys:
582: +     -snes_ms_type - type of multi-stage smoother
583: -     -snes_ms_damping - damping for multi-stage method

585:       Notes:
586:       These multistage methods are explicit Runge-Kutta methods that are often used as smoothers for
587:       FAS multigrid for transport problems. In the linear case, these are equivalent to polynomial smoothers (such as Chebyshev).

589:       Multi-stage smoothers should usually be preconditioned by point-block Jacobi to ensure proper scaling and to normalize the wave speeds.

591:       The methods are specified in low storage form (Ketcheson 2010). New methods can be registered with `SNESMSRegister()`.

593:       References:
594: +     * - Ketcheson (2010) Runge Kutta methods with minimum storage implementations (https://doi.org/10.1016/j.jcp.2009.11.006).
595: .     * - 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).
596: .     * - Pierce and Giles (1997) Preconditioned multigrid methods for compressible flow calculations on stretched meshes (https://doi.org/10.1006/jcph.1997.5772).
597: -     * - 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).

599:       Level: advanced

601: .seealso: `SNESCreate()`, `SNES`, `SNESSetType()`, `SNESMS`, `SNESFAS`, `KSPCHEBYSHEV`

603: M*/
604: PETSC_EXTERN PetscErrorCode SNESCreate_MS(SNES snes)
605: {
606:   SNES_MS *ms;

608:   PetscFunctionBegin;
609:   PetscCall(SNESMSInitializePackage());

611:   snes->ops->setup          = SNESSetUp_MS;
612:   snes->ops->solve          = SNESSolve_MS;
613:   snes->ops->destroy        = SNESDestroy_MS;
614:   snes->ops->setfromoptions = SNESSetFromOptions_MS;
615:   snes->ops->view           = SNESView_MS;
616:   snes->ops->reset          = SNESReset_MS;

618:   snes->usesnpc = PETSC_FALSE;
619:   snes->usesksp = PETSC_TRUE;

621:   snes->alwayscomputesfinalresidual = PETSC_FALSE;

623:   PetscCall(PetscNew(&ms));
624:   snes->data  = (void *)ms;
625:   ms->damping = 0.9;

627:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetType_C", SNESMSGetType_MS));
628:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetType_C", SNESMSSetType_MS));
629:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSGetDamping_C", SNESMSGetDamping_MS));
630:   PetscCall(PetscObjectComposeFunction((PetscObject)snes, "SNESMSSetDamping_C", SNESMSSetDamping_MS));

632:   PetscCall(SNESMSSetType(snes, SNESMSDefault));
633:   PetscFunctionReturn(PETSC_SUCCESS);
634: }