Actual source code: rosw.c

  1: /*
  2:   Code for timestepping with Rosenbrock W methods

  4:   Notes:
  5:   The general system is written as

  7:   F(t,U,Udot) = G(t,U)

  9:   where F represents the stiff part of the physics and G represents the non-stiff part.
 10:   This method is designed to be linearly implicit on F and can use an approximate and lagged Jacobian.

 12: */
 13: #include <petsc/private/tsimpl.h>
 14: #include <petscdm.h>

 16: #include <petsc/private/kernels/blockinvert.h>

 18: static TSRosWType TSRosWDefault = TSROSWRA34PW2;
 19: static PetscBool  TSRosWRegisterAllCalled;
 20: static PetscBool  TSRosWPackageInitialized;

 22: typedef struct _RosWTableau *RosWTableau;
 23: struct _RosWTableau {
 24:   char      *name;
 25:   PetscInt   order;             /* Classical approximation order of the method */
 26:   PetscInt   s;                 /* Number of stages */
 27:   PetscInt   pinterp;           /* Interpolation order */
 28:   PetscReal *A;                 /* Propagation table, strictly lower triangular */
 29:   PetscReal *Gamma;             /* Stage table, lower triangular with nonzero diagonal */
 30:   PetscBool *GammaZeroDiag;     /* Diagonal entries that are zero in stage table Gamma, vector indicating explicit statages */
 31:   PetscReal *GammaExplicitCorr; /* Coefficients for correction terms needed for explicit stages in transformed variables*/
 32:   PetscReal *b;                 /* Step completion table */
 33:   PetscReal *bembed;            /* Step completion table for embedded method of order one less */
 34:   PetscReal *ASum;              /* Row sum of A */
 35:   PetscReal *GammaSum;          /* Row sum of Gamma, only needed for non-autonomous systems */
 36:   PetscReal *At;                /* Propagation table in transformed variables */
 37:   PetscReal *bt;                /* Step completion table in transformed variables */
 38:   PetscReal *bembedt;           /* Step completion table of order one less in transformed variables */
 39:   PetscReal *GammaInv;          /* Inverse of Gamma, used for transformed variables */
 40:   PetscReal  ccfl;              /* Placeholder for CFL coefficient relative to forward Euler */
 41:   PetscReal *binterpt;          /* Dense output formula */
 42: };
 43: typedef struct _RosWTableauLink *RosWTableauLink;
 44: struct _RosWTableauLink {
 45:   struct _RosWTableau tab;
 46:   RosWTableauLink     next;
 47: };
 48: static RosWTableauLink RosWTableauList;

 50: typedef struct {
 51:   RosWTableau  tableau;
 52:   Vec         *Y;            /* States computed during the step, used to complete the step */
 53:   Vec          Ydot;         /* Work vector holding Ydot during residual evaluation */
 54:   Vec          Ystage;       /* Work vector for the state value at each stage */
 55:   Vec          Zdot;         /* Ydot = Zdot + shift*Y */
 56:   Vec          Zstage;       /* Y = Zstage + Y */
 57:   Vec          vec_sol_prev; /* Solution from the previous step (used for interpolation and rollback)*/
 58:   PetscScalar *work;         /* Scalar work space of length number of stages, used to prepare VecMAXPY() */
 59:   PetscReal    scoeff;       /* shift = scoeff/dt */
 60:   PetscReal    stage_time;
 61:   PetscReal    stage_explicit;     /* Flag indicates that the current stage is explicit */
 62:   PetscBool    recompute_jacobian; /* Recompute the Jacobian at each stage, default is to freeze the Jacobian at the start of each step */
 63:   TSStepStatus status;
 64: } TS_RosW;

 66: /*MC
 67:      TSROSWTHETA1 - One stage first order L-stable Rosenbrock-W scheme (aka theta method).

 69:      Only an approximate Jacobian is needed.

 71:      Level: intermediate

 73: .seealso: [](ch_ts), `TSROSW`
 74: M*/

 76: /*MC
 77:      TSROSWTHETA2 - One stage second order A-stable Rosenbrock-W scheme (aka theta method).

 79:      Only an approximate Jacobian is needed.

 81:      Level: intermediate

 83: .seealso: [](ch_ts), `TSROSW`
 84: M*/

 86: /*MC
 87:      TSROSW2M - Two stage second order L-stable Rosenbrock-W scheme.

 89:      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2P.

 91:      Level: intermediate

 93: .seealso: [](ch_ts), `TSROSW`
 94: M*/

 96: /*MC
 97:      TSROSW2P - Two stage second order L-stable Rosenbrock-W scheme.

 99:      Only an approximate Jacobian is needed. By default, it is only recomputed once per step. This method is a reflection of TSROSW2M.

101:      Level: intermediate

103: .seealso: [](ch_ts), `TSROSW`
104: M*/

106: /*MC
107:      TSROSWRA3PW - Three stage third order Rosenbrock-W scheme for PDAE of index 1.

109:      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.

111:      This is strongly A-stable with R(infty) = 0.73. The embedded method of order 2 is strongly A-stable with R(infty) = 0.73.

113:      Level: intermediate

115:      References:
116: .  * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.

118: .seealso: [](ch_ts), `TSROSW`
119: M*/

121: /*MC
122:      TSROSWRA34PW2 - Four stage third order L-stable Rosenbrock-W scheme for PDAE of index 1.

124:      Only an approximate Jacobian is needed. By default, it is only recomputed once per step.

126:      This is strongly A-stable with R(infty) = 0. The embedded method of order 2 is strongly A-stable with R(infty) = 0.48.

128:      Level: intermediate

130:      References:
131: .  * - Rang and Angermann, New Rosenbrock W methods of order 3 for partial differential algebraic equations of index 1, 2005.

133: .seealso: [](ch_ts), `TSROSW`
134: M*/

136: /*MC
137:      TSROSWRODAS3 - Four stage third order L-stable Rosenbrock scheme

139:      By default, the Jacobian is only recomputed once per step.

141:      Both the third order and embedded second order methods are stiffly accurate and L-stable.

143:      Level: intermediate

145:      References:
146: .  * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.

148: .seealso: [](ch_ts), `TSROSW`, `TSROSWSANDU3`
149: M*/

151: /*MC
152:      TSROSWSANDU3 - Three stage third order L-stable Rosenbrock scheme

154:      By default, the Jacobian is only recomputed once per step.

156:      The third order method is L-stable, but not stiffly accurate.
157:      The second order embedded method is strongly A-stable with R(infty) = 0.5.
158:      The internal stages are L-stable.
159:      This method is called ROS3 in the paper.

161:      Level: intermediate

163:      References:
164: .  * - Sandu et al, Benchmarking stiff ODE solvers for atmospheric chemistry problems II, Rosenbrock solvers, 1997.

166: .seealso: [](ch_ts), `TSROSW`, `TSROSWRODAS3`
167: M*/

169: /*MC
170:      TSROSWASSP3P3S1C - A-stable Rosenbrock-W method with SSP explicit part, third order, three stages

172:      By default, the Jacobian is only recomputed once per step.

174:      A-stable SPP explicit order 3, 3 stages, CFL 1 (eff = 1/3)

176:      Level: intermediate

178:      References:
179: . * - Emil Constantinescu

181: .seealso: [](ch_ts), `TSROSW`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `SSP`
182: M*/

184: /*MC
185:      TSROSWLASSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages

187:      By default, the Jacobian is only recomputed once per step.

189:      L-stable (A-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)

191:      Level: intermediate

193:      References:
194: . * - Emil Constantinescu

196: .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLLSSP3P4S2C`, `TSSSP`
197: M*/

199: /*MC
200:      TSROSWLLSSP3P4S2C - L-stable Rosenbrock-W method with SSP explicit part, third order, four stages

202:      By default, the Jacobian is only recomputed once per step.

204:      L-stable (L-stable embedded) SPP explicit order 3, 4 stages, CFL 2 (eff = 1/2)

206:      Level: intermediate

208:      References:
209: . * - Emil Constantinescu

211: .seealso: [](ch_ts), `TSROSW`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSSSP`
212: M*/

214: /*MC
215:      TSROSWGRK4T - four stage, fourth order Rosenbrock (not W) method from Kaps and Rentrop

217:      By default, the Jacobian is only recomputed once per step.

219:      A(89.3 degrees)-stable, |R(infty)| = 0.454.

221:      This method does not provide a dense output formula.

223:      Level: intermediate

225:      References:
226: +   * -  Kaps and Rentrop, Generalized Runge Kutta methods of order four with stepsize control for stiff ordinary differential equations, 1979.
227: -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

229:      Hairer's code ros4.f

231: .seealso: [](ch_ts), `TSROSW`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`
232: M*/

234: /*MC
235:      TSROSWSHAMP4 - four stage, fourth order Rosenbrock (not W) method from Shampine

237:      By default, the Jacobian is only recomputed once per step.

239:      A-stable, |R(infty)| = 1/3.

241:      This method does not provide a dense output formula.

243:      Level: intermediate

245:      References:
246: +   * -  Shampine, Implementation of Rosenbrock methods, 1982.
247: -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

249:      Hairer's code ros4.f

251: .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWVELDD4`, `TSROSW4L`
252: M*/

254: /*MC
255:      TSROSWVELDD4 - four stage, fourth order Rosenbrock (not W) method from van Veldhuizen

257:      By default, the Jacobian is only recomputed once per step.

259:      A(89.5 degrees)-stable, |R(infty)| = 0.24.

261:      This method does not provide a dense output formula.

263:      Level: intermediate

265:      References:
266: +   * -  van Veldhuizen, D stability and Kaps Rentrop methods, 1984.
267: -   * -  Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

269:      Hairer's code ros4.f

271: .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L`
272: M*/

274: /*MC
275:      TSROSW4L - four stage, fourth order Rosenbrock (not W) method

277:      By default, the Jacobian is only recomputed once per step.

279:      A-stable and L-stable

281:      This method does not provide a dense output formula.

283:      Level: intermediate

285:      References:
286: .  * -   Hairer and Wanner, Solving Ordinary Differential Equations II, Section 4 Table 7.2.

288:      Hairer's code ros4.f

290: .seealso: [](ch_ts), `TSROSW`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSW4L`
291: M*/

293: /*@C
294:   TSRosWRegisterAll - Registers all of the Rosenbrock-W methods in `TSROSW`

296:   Not Collective, but should be called by all processes which will need the schemes to be registered

298:   Level: advanced

300: .seealso: [](ch_ts), `TSROSW`, `TSRosWRegisterDestroy()`
301: @*/
302: PetscErrorCode TSRosWRegisterAll(void)
303: {
304:   PetscFunctionBegin;
305:   if (TSRosWRegisterAllCalled) PetscFunctionReturn(PETSC_SUCCESS);
306:   TSRosWRegisterAllCalled = PETSC_TRUE;

308:   {
309:     const PetscReal A        = 0;
310:     const PetscReal Gamma    = 1;
311:     const PetscReal b        = 1;
312:     const PetscReal binterpt = 1;

314:     PetscCall(TSRosWRegister(TSROSWTHETA1, 1, 1, &A, &Gamma, &b, NULL, 1, &binterpt));
315:   }

317:   {
318:     const PetscReal A        = 0;
319:     const PetscReal Gamma    = 0.5;
320:     const PetscReal b        = 1;
321:     const PetscReal binterpt = 1;

323:     PetscCall(TSRosWRegister(TSROSWTHETA2, 2, 1, &A, &Gamma, &b, NULL, 1, &binterpt));
324:   }

326:   {
327:     /*const PetscReal g = 1. + 1./PetscSqrtReal(2.0);   Direct evaluation: 1.707106781186547524401. Used for setting up arrays of values known at compile time below. */
328:     const PetscReal A[2][2] = {
329:       {0,  0},
330:       {1., 0}
331:     };
332:     const PetscReal Gamma[2][2] = {
333:       {1.707106781186547524401,       0                      },
334:       {-2. * 1.707106781186547524401, 1.707106781186547524401}
335:     };
336:     const PetscReal b[2]  = {0.5, 0.5};
337:     const PetscReal b1[2] = {1.0, 0.0};
338:     PetscReal       binterpt[2][2];
339:     binterpt[0][0] = 1.707106781186547524401 - 1.0;
340:     binterpt[1][0] = 2.0 - 1.707106781186547524401;
341:     binterpt[0][1] = 1.707106781186547524401 - 1.5;
342:     binterpt[1][1] = 1.5 - 1.707106781186547524401;

344:     PetscCall(TSRosWRegister(TSROSW2P, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
345:   }
346:   {
347:     /*const PetscReal g = 1. - 1./PetscSqrtReal(2.0);   Direct evaluation: 0.2928932188134524755992. Used for setting up arrays of values known at compile time below. */
348:     const PetscReal A[2][2] = {
349:       {0,  0},
350:       {1., 0}
351:     };
352:     const PetscReal Gamma[2][2] = {
353:       {0.2928932188134524755992,       0                       },
354:       {-2. * 0.2928932188134524755992, 0.2928932188134524755992}
355:     };
356:     const PetscReal b[2]  = {0.5, 0.5};
357:     const PetscReal b1[2] = {1.0, 0.0};
358:     PetscReal       binterpt[2][2];
359:     binterpt[0][0] = 0.2928932188134524755992 - 1.0;
360:     binterpt[1][0] = 2.0 - 0.2928932188134524755992;
361:     binterpt[0][1] = 0.2928932188134524755992 - 1.5;
362:     binterpt[1][1] = 1.5 - 0.2928932188134524755992;

364:     PetscCall(TSRosWRegister(TSROSW2M, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
365:   }
366:   {
367:     /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */
368:     PetscReal       binterpt[3][2];
369:     const PetscReal A[3][3] = {
370:       {0,                      0, 0},
371:       {1.5773502691896257e+00, 0, 0},
372:       {0.5,                    0, 0}
373:     };
374:     const PetscReal Gamma[3][3] = {
375:       {7.8867513459481287e-01,  0,                       0                     },
376:       {-1.5773502691896257e+00, 7.8867513459481287e-01,  0                     },
377:       {-6.7075317547305480e-01, -1.7075317547305482e-01, 7.8867513459481287e-01}
378:     };
379:     const PetscReal b[3]  = {1.0566243270259355e-01, 4.9038105676657971e-02, 8.4529946162074843e-01};
380:     const PetscReal b2[3] = {-1.7863279495408180e-01, 1. / 3., 8.4529946162074843e-01};

382:     binterpt[0][0] = -0.8094010767585034;
383:     binterpt[1][0] = -0.5;
384:     binterpt[2][0] = 2.3094010767585034;
385:     binterpt[0][1] = 0.9641016151377548;
386:     binterpt[1][1] = 0.5;
387:     binterpt[2][1] = -1.4641016151377548;

389:     PetscCall(TSRosWRegister(TSROSWRA3PW, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
390:   }
391:   {
392:     PetscReal binterpt[4][3];
393:     /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */
394:     const PetscReal A[4][4] = {
395:       {0,                      0,                       0,  0},
396:       {8.7173304301691801e-01, 0,                       0,  0},
397:       {8.4457060015369423e-01, -1.1299064236484185e-01, 0,  0},
398:       {0,                      0,                       1., 0}
399:     };
400:     const PetscReal Gamma[4][4] = {
401:       {4.3586652150845900e-01,  0,                       0,                      0                     },
402:       {-8.7173304301691801e-01, 4.3586652150845900e-01,  0,                      0                     },
403:       {-9.0338057013044082e-01, 5.4180672388095326e-02,  4.3586652150845900e-01, 0                     },
404:       {2.4212380706095346e-01,  -1.2232505839045147e+00, 5.4526025533510214e-01, 4.3586652150845900e-01}
405:     };
406:     const PetscReal b[4]  = {2.4212380706095346e-01, -1.2232505839045147e+00, 1.5452602553351020e+00, 4.3586652150845900e-01};
407:     const PetscReal b2[4] = {3.7810903145819369e-01, -9.6042292212423178e-02, 5.0000000000000000e-01, 2.1793326075422950e-01};

409:     binterpt[0][0] = 1.0564298455794094;
410:     binterpt[1][0] = 2.296429974281067;
411:     binterpt[2][0] = -1.307599564525376;
412:     binterpt[3][0] = -1.045260255335102;
413:     binterpt[0][1] = -1.3864882699759573;
414:     binterpt[1][1] = -8.262611700275677;
415:     binterpt[2][1] = 7.250979895056055;
416:     binterpt[3][1] = 2.398120075195581;
417:     binterpt[0][2] = 0.5721822314575016;
418:     binterpt[1][2] = 4.742931142090097;
419:     binterpt[2][2] = -4.398120075195578;
420:     binterpt[3][2] = -0.9169932983520199;

422:     PetscCall(TSRosWRegister(TSROSWRA34PW2, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
423:   }
424:   {
425:     /* const PetscReal g = 0.5;       Directly written in-place below */
426:     const PetscReal A[4][4] = {
427:       {0,    0,     0,   0},
428:       {0,    0,     0,   0},
429:       {1.,   0,     0,   0},
430:       {0.75, -0.25, 0.5, 0}
431:     };
432:     const PetscReal Gamma[4][4] = {
433:       {0.5,     0,       0,       0  },
434:       {1.,      0.5,     0,       0  },
435:       {-0.25,   -0.25,   0.5,     0  },
436:       {1. / 12, 1. / 12, -2. / 3, 0.5}
437:     };
438:     const PetscReal b[4]  = {5. / 6, -1. / 6, -1. / 6, 0.5};
439:     const PetscReal b2[4] = {0.75, -0.25, 0.5, 0};

441:     PetscCall(TSRosWRegister(TSROSWRODAS3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 0, NULL));
442:   }
443:   {
444:     /*const PetscReal g = 0.43586652150845899941601945119356;       Directly written in-place below */
445:     const PetscReal A[3][3] = {
446:       {0,                                  0, 0},
447:       {0.43586652150845899941601945119356, 0, 0},
448:       {0.43586652150845899941601945119356, 0, 0}
449:     };
450:     const PetscReal Gamma[3][3] = {
451:       {0.43586652150845899941601945119356,  0,                                  0                                 },
452:       {-0.19294655696029095575009695436041, 0.43586652150845899941601945119356, 0                                 },
453:       {0,                                   1.74927148125794685173529749738960, 0.43586652150845899941601945119356}
454:     };
455:     const PetscReal b[3]  = {-0.75457412385404315829818998646589, 1.94100407061964420292840123379419, -0.18642994676560104463021124732829};
456:     const PetscReal b2[3] = {-1.53358745784149585370766523913002, 2.81745131148625772213931745457622, -0.28386385364476186843165221544619};

458:     PetscReal binterpt[3][2];
459:     binterpt[0][0] = 3.793692883777660870425141387941;
460:     binterpt[1][0] = -2.918692883777660870425141387941;
461:     binterpt[2][0] = 0.125;
462:     binterpt[0][1] = -0.725741064379812106687651020584;
463:     binterpt[1][1] = 0.559074397713145440020984353917;
464:     binterpt[2][1] = 0.16666666666666666666666666666667;

466:     PetscCall(TSRosWRegister(TSROSWSANDU3, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
467:   }
468:   {
469:     /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
470:      * Direct evaluation: s3 = 1.732050807568877293527;
471:      *                     g = 0.7886751345948128822546;
472:      * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */
473:     const PetscReal A[3][3] = {
474:       {0,    0,    0},
475:       {1,    0,    0},
476:       {0.25, 0.25, 0}
477:     };
478:     const PetscReal Gamma[3][3] = {
479:       {0,                                       0,                                      0                       },
480:       {(-3.0 - 1.732050807568877293527) / 6.0,  0.7886751345948128822546,               0                       },
481:       {(-3.0 - 1.732050807568877293527) / 24.0, (-3.0 - 1.732050807568877293527) / 8.0, 0.7886751345948128822546}
482:     };
483:     const PetscReal b[3]  = {1. / 6., 1. / 6., 2. / 3.};
484:     const PetscReal b2[3] = {1. / 4., 1. / 4., 1. / 2.};
485:     PetscReal       binterpt[3][2];

487:     binterpt[0][0] = 0.089316397477040902157517886164709;
488:     binterpt[1][0] = -0.91068360252295909784248211383529;
489:     binterpt[2][0] = 1.8213672050459181956849642276706;
490:     binterpt[0][1] = 0.077350269189625764509148780501957;
491:     binterpt[1][1] = 1.077350269189625764509148780502;
492:     binterpt[2][1] = -1.1547005383792515290182975610039;

494:     PetscCall(TSRosWRegister(TSROSWASSP3P3S1C, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
495:   }

497:   {
498:     const PetscReal A[4][4] = {
499:       {0,       0,       0,       0},
500:       {1. / 2., 0,       0,       0},
501:       {1. / 2., 1. / 2., 0,       0},
502:       {1. / 6., 1. / 6., 1. / 6., 0}
503:     };
504:     const PetscReal Gamma[4][4] = {
505:       {1. / 2., 0,        0,       0},
506:       {0.0,     1. / 4.,  0,       0},
507:       {-2.,     -2. / 3., 2. / 3., 0},
508:       {1. / 2., 5. / 36., -2. / 9, 0}
509:     };
510:     const PetscReal b[4]  = {1. / 6., 1. / 6., 1. / 6., 1. / 2.};
511:     const PetscReal b2[4] = {1. / 8., 3. / 4., 1. / 8., 0};
512:     PetscReal       binterpt[4][3];

514:     binterpt[0][0] = 6.25;
515:     binterpt[1][0] = -30.25;
516:     binterpt[2][0] = 1.75;
517:     binterpt[3][0] = 23.25;
518:     binterpt[0][1] = -9.75;
519:     binterpt[1][1] = 58.75;
520:     binterpt[2][1] = -3.25;
521:     binterpt[3][1] = -45.75;
522:     binterpt[0][2] = 3.6666666666666666666666666666667;
523:     binterpt[1][2] = -28.333333333333333333333333333333;
524:     binterpt[2][2] = 1.6666666666666666666666666666667;
525:     binterpt[3][2] = 23.;

527:     PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
528:   }

530:   {
531:     const PetscReal A[4][4] = {
532:       {0,       0,       0,       0},
533:       {1. / 2., 0,       0,       0},
534:       {1. / 2., 1. / 2., 0,       0},
535:       {1. / 6., 1. / 6., 1. / 6., 0}
536:     };
537:     const PetscReal Gamma[4][4] = {
538:       {1. / 2.,  0,          0,        0},
539:       {0.0,      3. / 4.,    0,        0},
540:       {-2. / 3., -23. / 9.,  2. / 9.,  0},
541:       {1. / 18., 65. / 108., -2. / 27, 0}
542:     };
543:     const PetscReal b[4]  = {1. / 6., 1. / 6., 1. / 6., 1. / 2.};
544:     const PetscReal b2[4] = {3. / 16., 10. / 16., 3. / 16., 0};
545:     PetscReal       binterpt[4][3];

547:     binterpt[0][0] = 1.6911764705882352941176470588235;
548:     binterpt[1][0] = 3.6813725490196078431372549019608;
549:     binterpt[2][0] = 0.23039215686274509803921568627451;
550:     binterpt[3][0] = -4.6029411764705882352941176470588;
551:     binterpt[0][1] = -0.95588235294117647058823529411765;
552:     binterpt[1][1] = -6.2401960784313725490196078431373;
553:     binterpt[2][1] = -0.31862745098039215686274509803922;
554:     binterpt[3][1] = 7.5147058823529411764705882352941;
555:     binterpt[0][2] = -0.56862745098039215686274509803922;
556:     binterpt[1][2] = 2.7254901960784313725490196078431;
557:     binterpt[2][2] = 0.25490196078431372549019607843137;
558:     binterpt[3][2] = -2.4117647058823529411764705882353;

560:     PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
561:   }

563:   {
564:     PetscReal A[4][4], Gamma[4][4], b[4], b2[4];
565:     PetscReal binterpt[4][3];

567:     Gamma[0][0] = 0.4358665215084589994160194475295062513822671686978816;
568:     Gamma[0][1] = 0;
569:     Gamma[0][2] = 0;
570:     Gamma[0][3] = 0;
571:     Gamma[1][0] = -1.997527830934941248426324674704153457289527280554476;
572:     Gamma[1][1] = 0.4358665215084589994160194475295062513822671686978816;
573:     Gamma[1][2] = 0;
574:     Gamma[1][3] = 0;
575:     Gamma[2][0] = -1.007948511795029620852002345345404191008352770119903;
576:     Gamma[2][1] = -0.004648958462629345562774289390054679806993396798458131;
577:     Gamma[2][2] = 0.4358665215084589994160194475295062513822671686978816;
578:     Gamma[2][3] = 0;
579:     Gamma[3][0] = -0.6685429734233467180451604600279552604364311322650783;
580:     Gamma[3][1] = 0.6056625986449338476089525334450053439525178740492984;
581:     Gamma[3][2] = -0.9717899277217721234705114616271378792182450260943198;
582:     Gamma[3][3] = 0;

584:     A[0][0] = 0;
585:     A[0][1] = 0;
586:     A[0][2] = 0;
587:     A[0][3] = 0;
588:     A[1][0] = 0.8717330430169179988320388950590125027645343373957631;
589:     A[1][1] = 0;
590:     A[1][2] = 0;
591:     A[1][3] = 0;
592:     A[2][0] = 0.5275890119763004115618079766722914408876108660811028;
593:     A[2][1] = 0.07241098802369958843819203208518599088698057726988732;
594:     A[2][2] = 0;
595:     A[2][3] = 0;
596:     A[3][0] = 0.3990960076760701320627260685975778145384666450351314;
597:     A[3][1] = -0.4375576546135194437228463747348862825846903771419953;
598:     A[3][2] = 1.038461646937449311660120300601880176655352737312713;
599:     A[3][3] = 0;

601:     b[0] = 0.1876410243467238251612921333138006734899663569186926;
602:     b[1] = -0.5952974735769549480478230473706443582188442040780541;
603:     b[2] = 0.9717899277217721234705114616271378792182450260943198;
604:     b[3] = 0.4358665215084589994160194475295062513822671686978816;

606:     b2[0] = 0.2147402862233891404862383521089097657790734483804460;
607:     b2[1] = -0.4851622638849390928209050538171743017757490232519684;
608:     b2[2] = 0.8687250025203875511662123688667549217531982787600080;
609:     b2[3] = 0.4016969751411624011684543450940068201770721128357014;

611:     binterpt[0][0] = 2.2565812720167954547104627844105;
612:     binterpt[1][0] = 1.349166413351089573796243820819;
613:     binterpt[2][0] = -2.4695174540533503758652847586647;
614:     binterpt[3][0] = -0.13623023131453465264142184656474;
615:     binterpt[0][1] = -3.0826699111559187902922463354557;
616:     binterpt[1][1] = -2.4689115685996042534544925650515;
617:     binterpt[2][1] = 5.7428279814696677152129332773553;
618:     binterpt[3][1] = -0.19124650171414467146619437684812;
619:     binterpt[0][2] = 1.0137296634858471607430756831148;
620:     binterpt[1][2] = 0.52444768167155973161042570784064;
621:     binterpt[2][2] = -2.3015205996945452158771370439586;
622:     binterpt[3][2] = 0.76334325453713832352363565300308;

624:     PetscCall(TSRosWRegister(TSROSWARK3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
625:   }
626:   PetscCall(TSRosWRegisterRos4(TSROSWGRK4T, 0.231, PETSC_DEFAULT, PETSC_DEFAULT, 0, -0.1282612945269037e+01));
627:   PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4, 0.5, PETSC_DEFAULT, PETSC_DEFAULT, 0, 125. / 108.));
628:   PetscCall(TSRosWRegisterRos4(TSROSWVELDD4, 0.22570811482256823492, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.355958941201148));
629:   PetscCall(TSRosWRegisterRos4(TSROSW4L, 0.57282, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.093502252409163));
630:   PetscFunctionReturn(PETSC_SUCCESS);
631: }

633: /*@C
634:   TSRosWRegisterDestroy - Frees the list of schemes that were registered by `TSRosWRegister()`.

636:   Not Collective

638:   Level: advanced

640: .seealso: [](ch_ts), `TSRosWRegister()`, `TSRosWRegisterAll()`
641: @*/
642: PetscErrorCode TSRosWRegisterDestroy(void)
643: {
644:   RosWTableauLink link;

646:   PetscFunctionBegin;
647:   while ((link = RosWTableauList)) {
648:     RosWTableau t   = &link->tab;
649:     RosWTableauList = link->next;
650:     PetscCall(PetscFree5(t->A, t->Gamma, t->b, t->ASum, t->GammaSum));
651:     PetscCall(PetscFree5(t->At, t->bt, t->GammaInv, t->GammaZeroDiag, t->GammaExplicitCorr));
652:     PetscCall(PetscFree2(t->bembed, t->bembedt));
653:     PetscCall(PetscFree(t->binterpt));
654:     PetscCall(PetscFree(t->name));
655:     PetscCall(PetscFree(link));
656:   }
657:   TSRosWRegisterAllCalled = PETSC_FALSE;
658:   PetscFunctionReturn(PETSC_SUCCESS);
659: }

661: /*@C
662:   TSRosWInitializePackage - This function initializes everything in the `TSROSW` package. It is called
663:   from `TSInitializePackage()`.

665:   Level: developer

667: .seealso: [](ch_ts), `TSROSW`, `PetscInitialize()`, `TSRosWFinalizePackage()`
668: @*/
669: PetscErrorCode TSRosWInitializePackage(void)
670: {
671:   PetscFunctionBegin;
672:   if (TSRosWPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
673:   TSRosWPackageInitialized = PETSC_TRUE;
674:   PetscCall(TSRosWRegisterAll());
675:   PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage));
676:   PetscFunctionReturn(PETSC_SUCCESS);
677: }

679: /*@C
680:   TSRosWFinalizePackage - This function destroys everything in the `TSROSW` package. It is
681:   called from `PetscFinalize()`.

683:   Level: developer

685: .seealso: [](ch_ts), `TSROSW`, `PetscFinalize()`, `TSRosWInitializePackage()`
686: @*/
687: PetscErrorCode TSRosWFinalizePackage(void)
688: {
689:   PetscFunctionBegin;
690:   TSRosWPackageInitialized = PETSC_FALSE;
691:   PetscCall(TSRosWRegisterDestroy());
692:   PetscFunctionReturn(PETSC_SUCCESS);
693: }

695: /*@C
696:   TSRosWRegister - register a `TSROSW`, Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation

698:   Not Collective, but the same schemes should be registered on all processes on which they will be used

700:   Input Parameters:
701: + name     - identifier for method
702: . order    - approximation order of method
703: . s        - number of stages, this is the dimension of the matrices below
704: . A        - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
705: . Gamma    - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
706: . b        - Step completion table (dimension s)
707: . bembed   - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
708: . pinterp  - Order of the interpolation scheme, equal to the number of columns of binterpt
709: - binterpt - Coefficients of the interpolation formula (dimension s*pinterp)

711:   Level: advanced

713:   Note:
714:   Several Rosenbrock W methods are provided, this function is only needed to create new methods.

716: .seealso: [](ch_ts), `TSROSW`
717: @*/
718: PetscErrorCode TSRosWRegister(TSRosWType name, PetscInt order, PetscInt s, const PetscReal A[], const PetscReal Gamma[], const PetscReal b[], const PetscReal bembed[], PetscInt pinterp, const PetscReal binterpt[])
719: {
720:   RosWTableauLink link;
721:   RosWTableau     t;
722:   PetscInt        i, j, k;
723:   PetscScalar    *GammaInv;

725:   PetscFunctionBegin;
726:   PetscAssertPointer(name, 1);
727:   PetscAssertPointer(A, 4);
728:   PetscAssertPointer(Gamma, 5);
729:   PetscAssertPointer(b, 6);
730:   if (bembed) PetscAssertPointer(bembed, 7);

732:   PetscCall(TSRosWInitializePackage());
733:   PetscCall(PetscNew(&link));
734:   t = &link->tab;
735:   PetscCall(PetscStrallocpy(name, &t->name));
736:   t->order = order;
737:   t->s     = s;
738:   PetscCall(PetscMalloc5(s * s, &t->A, s * s, &t->Gamma, s, &t->b, s, &t->ASum, s, &t->GammaSum));
739:   PetscCall(PetscMalloc5(s * s, &t->At, s, &t->bt, s * s, &t->GammaInv, s, &t->GammaZeroDiag, s * s, &t->GammaExplicitCorr));
740:   PetscCall(PetscArraycpy(t->A, A, s * s));
741:   PetscCall(PetscArraycpy(t->Gamma, Gamma, s * s));
742:   PetscCall(PetscArraycpy(t->GammaExplicitCorr, Gamma, s * s));
743:   PetscCall(PetscArraycpy(t->b, b, s));
744:   if (bembed) {
745:     PetscCall(PetscMalloc2(s, &t->bembed, s, &t->bembedt));
746:     PetscCall(PetscArraycpy(t->bembed, bembed, s));
747:   }
748:   for (i = 0; i < s; i++) {
749:     t->ASum[i]     = 0;
750:     t->GammaSum[i] = 0;
751:     for (j = 0; j < s; j++) {
752:       t->ASum[i] += A[i * s + j];
753:       t->GammaSum[i] += Gamma[i * s + j];
754:     }
755:   }
756:   PetscCall(PetscMalloc1(s * s, &GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */
757:   for (i = 0; i < s * s; i++) GammaInv[i] = Gamma[i];
758:   for (i = 0; i < s; i++) {
759:     if (Gamma[i * s + i] == 0.0) {
760:       GammaInv[i * s + i] = 1.0;
761:       t->GammaZeroDiag[i] = PETSC_TRUE;
762:     } else {
763:       t->GammaZeroDiag[i] = PETSC_FALSE;
764:     }
765:   }

767:   switch (s) {
768:   case 1:
769:     GammaInv[0] = 1. / GammaInv[0];
770:     break;
771:   case 2:
772:     PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv, 0, PETSC_FALSE, NULL));
773:     break;
774:   case 3:
775:     PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv, 0, PETSC_FALSE, NULL));
776:     break;
777:   case 4:
778:     PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv, 0, PETSC_FALSE, NULL));
779:     break;
780:   case 5: {
781:     PetscInt  ipvt5[5];
782:     MatScalar work5[5 * 5];
783:     PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv, ipvt5, work5, 0, PETSC_FALSE, NULL));
784:     break;
785:   }
786:   case 6:
787:     PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv, 0, PETSC_FALSE, NULL));
788:     break;
789:   case 7:
790:     PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv, 0, PETSC_FALSE, NULL));
791:     break;
792:   default:
793:     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " stages", s);
794:   }
795:   for (i = 0; i < s * s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
796:   PetscCall(PetscFree(GammaInv));

798:   for (i = 0; i < s; i++) {
799:     for (k = 0; k < i + 1; k++) {
800:       t->GammaExplicitCorr[i * s + k] = (t->GammaExplicitCorr[i * s + k]) * (t->GammaInv[k * s + k]);
801:       for (j = k + 1; j < i + 1; j++) t->GammaExplicitCorr[i * s + k] += (t->GammaExplicitCorr[i * s + j]) * (t->GammaInv[j * s + k]);
802:     }
803:   }

805:   for (i = 0; i < s; i++) {
806:     for (j = 0; j < s; j++) {
807:       t->At[i * s + j] = 0;
808:       for (k = 0; k < s; k++) t->At[i * s + j] += t->A[i * s + k] * t->GammaInv[k * s + j];
809:     }
810:     t->bt[i] = 0;
811:     for (j = 0; j < s; j++) t->bt[i] += t->b[j] * t->GammaInv[j * s + i];
812:     if (bembed) {
813:       t->bembedt[i] = 0;
814:       for (j = 0; j < s; j++) t->bembedt[i] += t->bembed[j] * t->GammaInv[j * s + i];
815:     }
816:   }
817:   t->ccfl = 1.0; /* Fix this */

819:   t->pinterp = pinterp;
820:   PetscCall(PetscMalloc1(s * pinterp, &t->binterpt));
821:   PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp));
822:   link->next      = RosWTableauList;
823:   RosWTableauList = link;
824:   PetscFunctionReturn(PETSC_SUCCESS);
825: }

827: /*@C
828:   TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices

830:   Not Collective, but the same schemes should be registered on all processes on which they will be used

832:   Input Parameters:
833: + name  - identifier for method
834: . gamma - leading coefficient (diagonal entry)
835: . a2    - design parameter, see Table 7.2 of Hairer&Wanner
836: . a3    - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
837: . b3    - design parameter, see Table 7.2 of Hairer&Wanner
838: - e4    - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer

840:   Level: developer

842:   Notes:
843:   This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2.
844:   It is used here to implement several methods from the book and can be used to experiment with new methods.
845:   It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions.

847: .seealso: [](ch_ts), `TSRosW`, `TSRosWRegister()`
848: @*/
849: PetscErrorCode TSRosWRegisterRos4(TSRosWType name, PetscReal gamma, PetscReal a2, PetscReal a3, PetscReal b3, PetscReal e4)
850: {
851:   /* Declare numeric constants so they can be quad precision without being truncated at double */
852:   const PetscReal one = 1, two = 2, three = 3, four = 4, five = 5, six = 6, eight = 8, twelve = 12, twenty = 20, twentyfour = 24, p32 = one / six - gamma + gamma * gamma, p42 = one / eight - gamma / three, p43 = one / twelve - gamma / three, p44 = one / twentyfour - gamma / two + three / two * gamma * gamma - gamma * gamma * gamma, p56 = one / twenty - gamma / four;
853:   PetscReal   a4, a32, a42, a43, b1, b2, b4, beta2p, beta3p, beta4p, beta32, beta42, beta43, beta32beta2p, beta4jbetajp;
854:   PetscReal   A[4][4], Gamma[4][4], b[4], bm[4];
855:   PetscScalar M[3][3], rhs[3];

857:   PetscFunctionBegin;
858:   /* Step 1: choose Gamma (input) */
859:   /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */
860:   if (a3 == (PetscReal)PETSC_DEFAULT) a3 = (one / five - a2 / four) / (one / four - a2 / three); /* Eq 7.22 */
861:   a4 = a3;                                                                                       /* consequence of 7.20 */

863:   /* Solve order conditions 7.15a, 7.15c, 7.15e */
864:   M[0][0] = one;
865:   M[0][1] = one;
866:   M[0][2] = one; /* 7.15a */
867:   M[1][0] = 0.0;
868:   M[1][1] = a2 * a2;
869:   M[1][2] = a4 * a4; /* 7.15c */
870:   M[2][0] = 0.0;
871:   M[2][1] = a2 * a2 * a2;
872:   M[2][2] = a4 * a4 * a4; /* 7.15e */
873:   rhs[0]  = one - b3;
874:   rhs[1]  = one / three - a3 * a3 * b3;
875:   rhs[2]  = one / four - a3 * a3 * a3 * b3;
876:   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
877:   b1 = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
878:   b2 = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
879:   b4 = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);

881:   /* Step 3 */
882:   beta43       = (p56 - a2 * p43) / (b4 * a3 * a3 * (a3 - a2)); /* 7.21 */
883:   beta32beta2p = p44 / (b4 * beta43);                           /* 7.15h */
884:   beta4jbetajp = (p32 - b3 * beta32beta2p) / b4;
885:   M[0][0]      = b2;
886:   M[0][1]      = b3;
887:   M[0][2]      = b4;
888:   M[1][0]      = a4 * a4 * beta32beta2p - a3 * a3 * beta4jbetajp;
889:   M[1][1]      = a2 * a2 * beta4jbetajp;
890:   M[1][2]      = -a2 * a2 * beta32beta2p;
891:   M[2][0]      = b4 * beta43 * a3 * a3 - p43;
892:   M[2][1]      = -b4 * beta43 * a2 * a2;
893:   M[2][2]      = 0;
894:   rhs[0]       = one / two - gamma;
895:   rhs[1]       = 0;
896:   rhs[2]       = -a2 * a2 * p32;
897:   PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
898:   beta2p = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
899:   beta3p = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
900:   beta4p = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);

902:   /* Step 4: back-substitute */
903:   beta32 = beta32beta2p / beta2p;
904:   beta42 = (beta4jbetajp - beta43 * beta3p) / beta2p;

906:   /* Step 5: 7.15f and 7.20, then 7.16 */
907:   a43 = 0;
908:   a32 = p42 / (b3 * a3 * beta2p + b4 * a4 * beta2p);
909:   a42 = a32;

911:   A[0][0]     = 0;
912:   A[0][1]     = 0;
913:   A[0][2]     = 0;
914:   A[0][3]     = 0;
915:   A[1][0]     = a2;
916:   A[1][1]     = 0;
917:   A[1][2]     = 0;
918:   A[1][3]     = 0;
919:   A[2][0]     = a3 - a32;
920:   A[2][1]     = a32;
921:   A[2][2]     = 0;
922:   A[2][3]     = 0;
923:   A[3][0]     = a4 - a43 - a42;
924:   A[3][1]     = a42;
925:   A[3][2]     = a43;
926:   A[3][3]     = 0;
927:   Gamma[0][0] = gamma;
928:   Gamma[0][1] = 0;
929:   Gamma[0][2] = 0;
930:   Gamma[0][3] = 0;
931:   Gamma[1][0] = beta2p - A[1][0];
932:   Gamma[1][1] = gamma;
933:   Gamma[1][2] = 0;
934:   Gamma[1][3] = 0;
935:   Gamma[2][0] = beta3p - beta32 - A[2][0];
936:   Gamma[2][1] = beta32 - A[2][1];
937:   Gamma[2][2] = gamma;
938:   Gamma[2][3] = 0;
939:   Gamma[3][0] = beta4p - beta42 - beta43 - A[3][0];
940:   Gamma[3][1] = beta42 - A[3][1];
941:   Gamma[3][2] = beta43 - A[3][2];
942:   Gamma[3][3] = gamma;
943:   b[0]        = b1;
944:   b[1]        = b2;
945:   b[2]        = b3;
946:   b[3]        = b4;

948:   /* Construct embedded formula using given e4. We are solving Equation 7.18. */
949:   bm[3] = b[3] - e4 * gamma;                                              /* using definition of E4 */
950:   bm[2] = (p32 - beta4jbetajp * bm[3]) / (beta32 * beta2p);               /* fourth row of 7.18 */
951:   bm[1] = (one / two - gamma - beta3p * bm[2] - beta4p * bm[3]) / beta2p; /* second row */
952:   bm[0] = one - bm[1] - bm[2] - bm[3];                                    /* first row */

954:   {
955:     const PetscReal misfit = a2 * a2 * bm[1] + a3 * a3 * bm[2] + a4 * a4 * bm[3] - one / three;
956:     PetscCheck(PetscAbs(misfit) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_SUP, "Assumptions violated, could not construct a third order embedded method");
957:   }
958:   PetscCall(TSRosWRegister(name, 4, 4, &A[0][0], &Gamma[0][0], b, bm, 0, NULL));
959:   PetscFunctionReturn(PETSC_SUCCESS);
960: }

962: /*
963:  The step completion formula is

965:  x1 = x0 + b^T Y

967:  where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
968:  updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write

970:  x1e = x0 + be^T Y
971:      = x1 - b^T Y + be^T Y
972:      = x1 + (be - b)^T Y

974:  so we can evaluate the method of different order even after the step has been optimistically completed.
975: */
976: static PetscErrorCode TSEvaluateStep_RosW(TS ts, PetscInt order, Vec U, PetscBool *done)
977: {
978:   TS_RosW     *ros = (TS_RosW *)ts->data;
979:   RosWTableau  tab = ros->tableau;
980:   PetscScalar *w   = ros->work;
981:   PetscInt     i;

983:   PetscFunctionBegin;
984:   if (order == tab->order) {
985:     if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
986:       PetscCall(VecCopy(ts->vec_sol, U));
987:       for (i = 0; i < tab->s; i++) w[i] = tab->bt[i];
988:       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
989:     } else PetscCall(VecCopy(ts->vec_sol, U));
990:     if (done) *done = PETSC_TRUE;
991:     PetscFunctionReturn(PETSC_SUCCESS);
992:   } else if (order == tab->order - 1) {
993:     if (!tab->bembedt) goto unavailable;
994:     if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
995:       PetscCall(VecCopy(ts->vec_sol, U));
996:       for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i];
997:       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
998:     } else { /* Use rollback-and-recomplete formula (bembedt - bt) */
999:       for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
1000:       PetscCall(VecCopy(ts->vec_sol, U));
1001:       PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
1002:     }
1003:     if (done) *done = PETSC_TRUE;
1004:     PetscFunctionReturn(PETSC_SUCCESS);
1005:   }
1006: unavailable:
1007:   if (done) *done = PETSC_FALSE;
1008:   else
1009:     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "Rosenbrock-W '%s' of order %" PetscInt_FMT " cannot evaluate step at order %" PetscInt_FMT ". Consider using -ts_adapt_type none or a different method that has an embedded estimate.", tab->name,
1010:             tab->order, order);
1011:   PetscFunctionReturn(PETSC_SUCCESS);
1012: }

1014: static PetscErrorCode TSRollBack_RosW(TS ts)
1015: {
1016:   TS_RosW *ros = (TS_RosW *)ts->data;

1018:   PetscFunctionBegin;
1019:   PetscCall(VecCopy(ros->vec_sol_prev, ts->vec_sol));
1020:   PetscFunctionReturn(PETSC_SUCCESS);
1021: }

1023: static PetscErrorCode TSStep_RosW(TS ts)
1024: {
1025:   TS_RosW         *ros = (TS_RosW *)ts->data;
1026:   RosWTableau      tab = ros->tableau;
1027:   const PetscInt   s   = tab->s;
1028:   const PetscReal *At = tab->At, *Gamma = tab->Gamma, *ASum = tab->ASum, *GammaInv = tab->GammaInv;
1029:   const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
1030:   const PetscBool *GammaZeroDiag     = tab->GammaZeroDiag;
1031:   PetscScalar     *w                 = ros->work;
1032:   Vec             *Y = ros->Y, Ydot = ros->Ydot, Zdot = ros->Zdot, Zstage = ros->Zstage;
1033:   SNES             snes;
1034:   TSAdapt          adapt;
1035:   PetscInt         i, j, its, lits;
1036:   PetscInt         rejections = 0;
1037:   PetscBool        stageok, accept = PETSC_TRUE;
1038:   PetscReal        next_time_step = ts->time_step;
1039:   PetscInt         lag;

1041:   PetscFunctionBegin;
1042:   if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, ros->vec_sol_prev));

1044:   ros->status = TS_STEP_INCOMPLETE;
1045:   while (!ts->reason && ros->status != TS_STEP_COMPLETE) {
1046:     const PetscReal h = ts->time_step;
1047:     for (i = 0; i < s; i++) {
1048:       ros->stage_time = ts->ptime + h * ASum[i];
1049:       PetscCall(TSPreStage(ts, ros->stage_time));
1050:       if (GammaZeroDiag[i]) {
1051:         ros->stage_explicit = PETSC_TRUE;
1052:         ros->scoeff         = 1.;
1053:       } else {
1054:         ros->stage_explicit = PETSC_FALSE;
1055:         ros->scoeff         = 1. / Gamma[i * s + i];
1056:       }

1058:       PetscCall(VecCopy(ts->vec_sol, Zstage));
1059:       for (j = 0; j < i; j++) w[j] = At[i * s + j];
1060:       PetscCall(VecMAXPY(Zstage, i, w, Y));

1062:       for (j = 0; j < i; j++) w[j] = 1. / h * GammaInv[i * s + j];
1063:       PetscCall(VecZeroEntries(Zdot));
1064:       PetscCall(VecMAXPY(Zdot, i, w, Y));

1066:       /* Initial guess taken from last stage */
1067:       PetscCall(VecZeroEntries(Y[i]));

1069:       if (!ros->stage_explicit) {
1070:         PetscCall(TSGetSNES(ts, &snes));
1071:         if (!ros->recompute_jacobian && !i) {
1072:           PetscCall(SNESGetLagJacobian(snes, &lag));
1073:           if (lag == 1) {                            /* use did not set a nontrivial lag, so lag over all stages */
1074:             PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */
1075:           }
1076:         }
1077:         PetscCall(SNESSolve(snes, NULL, Y[i]));
1078:         if (!ros->recompute_jacobian && i == s - 1 && lag == 1) { PetscCall(SNESSetLagJacobian(snes, lag)); /* Set lag back to 1 so we know user did not set it */ }
1079:         PetscCall(SNESGetIterationNumber(snes, &its));
1080:         PetscCall(SNESGetLinearSolveIterations(snes, &lits));
1081:         ts->snes_its += its;
1082:         ts->ksp_its += lits;
1083:       } else {
1084:         Mat J, Jp;
1085:         PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
1086:         PetscCall(TSComputeIFunction(ts, ros->stage_time, Zstage, Ydot, Y[i], PETSC_FALSE));
1087:         PetscCall(VecScale(Y[i], -1.0));
1088:         PetscCall(VecAXPY(Y[i], -1.0, Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/

1090:         PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
1091:         for (j = 0; j < i; j++) w[j] = GammaExplicitCorr[i * s + j];
1092:         PetscCall(VecMAXPY(Zstage, i, w, Y));

1094:         /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
1095:         PetscCall(TSGetIJacobian(ts, &J, &Jp, NULL, NULL));
1096:         PetscCall(TSComputeIJacobian(ts, ros->stage_time, ts->vec_sol, Ydot, 0, J, Jp, PETSC_FALSE));
1097:         PetscCall(MatMult(J, Zstage, Zdot));
1098:         PetscCall(VecAXPY(Y[i], -1.0, Zdot));
1099:         ts->ksp_its += 1;

1101:         PetscCall(VecScale(Y[i], h));
1102:       }
1103:       PetscCall(TSPostStage(ts, ros->stage_time, i, Y));
1104:       PetscCall(TSGetAdapt(ts, &adapt));
1105:       PetscCall(TSAdaptCheckStage(adapt, ts, ros->stage_time, Y[i], &stageok));
1106:       if (!stageok) goto reject_step;
1107:     }

1109:     ros->status = TS_STEP_INCOMPLETE;
1110:     PetscCall(TSEvaluateStep_RosW(ts, tab->order, ts->vec_sol, NULL));
1111:     ros->status = TS_STEP_PENDING;
1112:     PetscCall(TSGetAdapt(ts, &adapt));
1113:     PetscCall(TSAdaptCandidatesClear(adapt));
1114:     PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
1115:     PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
1116:     ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
1117:     if (!accept) { /* Roll back the current step */
1118:       PetscCall(TSRollBack_RosW(ts));
1119:       ts->time_step = next_time_step;
1120:       goto reject_step;
1121:     }

1123:     ts->ptime += ts->time_step;
1124:     ts->time_step = next_time_step;
1125:     break;

1127:   reject_step:
1128:     ts->reject++;
1129:     accept = PETSC_FALSE;
1130:     if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
1131:       ts->reason = TS_DIVERGED_STEP_REJECTED;
1132:       PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
1133:     }
1134:   }
1135:   PetscFunctionReturn(PETSC_SUCCESS);
1136: }

1138: static PetscErrorCode TSInterpolate_RosW(TS ts, PetscReal itime, Vec U)
1139: {
1140:   TS_RosW         *ros = (TS_RosW *)ts->data;
1141:   PetscInt         s = ros->tableau->s, pinterp = ros->tableau->pinterp, i, j;
1142:   PetscReal        h;
1143:   PetscReal        tt, t;
1144:   PetscScalar     *bt;
1145:   const PetscReal *Bt       = ros->tableau->binterpt;
1146:   const PetscReal *GammaInv = ros->tableau->GammaInv;
1147:   PetscScalar     *w        = ros->work;
1148:   Vec             *Y        = ros->Y;

1150:   PetscFunctionBegin;
1151:   PetscCheck(Bt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRosW %s does not have an interpolation formula", ros->tableau->name);

1153:   switch (ros->status) {
1154:   case TS_STEP_INCOMPLETE:
1155:   case TS_STEP_PENDING:
1156:     h = ts->time_step;
1157:     t = (itime - ts->ptime) / h;
1158:     break;
1159:   case TS_STEP_COMPLETE:
1160:     h = ts->ptime - ts->ptime_prev;
1161:     t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
1162:     break;
1163:   default:
1164:     SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
1165:   }
1166:   PetscCall(PetscMalloc1(s, &bt));
1167:   for (i = 0; i < s; i++) bt[i] = 0;
1168:   for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
1169:     for (i = 0; i < s; i++) bt[i] += Bt[i * pinterp + j] * tt;
1170:   }

1172:   /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1173:   /* U <- 0*/
1174:   PetscCall(VecZeroEntries(U));
1175:   /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
1176:   for (j = 0; j < s; j++) w[j] = 0;
1177:   for (j = 0; j < s; j++) {
1178:     for (i = j; i < s; i++) w[j] += bt[i] * GammaInv[i * s + j];
1179:   }
1180:   PetscCall(VecMAXPY(U, i, w, Y));
1181:   /* U <- y(t) + U */
1182:   PetscCall(VecAXPY(U, 1, ros->vec_sol_prev));

1184:   PetscCall(PetscFree(bt));
1185:   PetscFunctionReturn(PETSC_SUCCESS);
1186: }

1188: /*------------------------------------------------------------*/

1190: static PetscErrorCode TSRosWTableauReset(TS ts)
1191: {
1192:   TS_RosW    *ros = (TS_RosW *)ts->data;
1193:   RosWTableau tab = ros->tableau;

1195:   PetscFunctionBegin;
1196:   if (!tab) PetscFunctionReturn(PETSC_SUCCESS);
1197:   PetscCall(VecDestroyVecs(tab->s, &ros->Y));
1198:   PetscCall(PetscFree(ros->work));
1199:   PetscFunctionReturn(PETSC_SUCCESS);
1200: }

1202: static PetscErrorCode TSReset_RosW(TS ts)
1203: {
1204:   TS_RosW *ros = (TS_RosW *)ts->data;

1206:   PetscFunctionBegin;
1207:   PetscCall(TSRosWTableauReset(ts));
1208:   PetscCall(VecDestroy(&ros->Ydot));
1209:   PetscCall(VecDestroy(&ros->Ystage));
1210:   PetscCall(VecDestroy(&ros->Zdot));
1211:   PetscCall(VecDestroy(&ros->Zstage));
1212:   PetscCall(VecDestroy(&ros->vec_sol_prev));
1213:   PetscFunctionReturn(PETSC_SUCCESS);
1214: }

1216: static PetscErrorCode TSRosWGetVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1217: {
1218:   TS_RosW *rw = (TS_RosW *)ts->data;

1220:   PetscFunctionBegin;
1221:   if (Ydot) {
1222:     if (dm && dm != ts->dm) {
1223:       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1224:     } else *Ydot = rw->Ydot;
1225:   }
1226:   if (Zdot) {
1227:     if (dm && dm != ts->dm) {
1228:       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1229:     } else *Zdot = rw->Zdot;
1230:   }
1231:   if (Ystage) {
1232:     if (dm && dm != ts->dm) {
1233:       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1234:     } else *Ystage = rw->Ystage;
1235:   }
1236:   if (Zstage) {
1237:     if (dm && dm != ts->dm) {
1238:       PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1239:     } else *Zstage = rw->Zstage;
1240:   }
1241:   PetscFunctionReturn(PETSC_SUCCESS);
1242: }

1244: static PetscErrorCode TSRosWRestoreVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1245: {
1246:   PetscFunctionBegin;
1247:   if (Ydot) {
1248:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1249:   }
1250:   if (Zdot) {
1251:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1252:   }
1253:   if (Ystage) {
1254:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1255:   }
1256:   if (Zstage) {
1257:     if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1258:   }
1259:   PetscFunctionReturn(PETSC_SUCCESS);
1260: }

1262: static PetscErrorCode DMCoarsenHook_TSRosW(DM fine, DM coarse, void *ctx)
1263: {
1264:   PetscFunctionBegin;
1265:   PetscFunctionReturn(PETSC_SUCCESS);
1266: }

1268: static PetscErrorCode DMRestrictHook_TSRosW(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
1269: {
1270:   TS  ts = (TS)ctx;
1271:   Vec Ydot, Zdot, Ystage, Zstage;
1272:   Vec Ydotc, Zdotc, Ystagec, Zstagec;

1274:   PetscFunctionBegin;
1275:   PetscCall(TSRosWGetVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
1276:   PetscCall(TSRosWGetVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
1277:   PetscCall(MatRestrict(restrct, Ydot, Ydotc));
1278:   PetscCall(VecPointwiseMult(Ydotc, rscale, Ydotc));
1279:   PetscCall(MatRestrict(restrct, Ystage, Ystagec));
1280:   PetscCall(VecPointwiseMult(Ystagec, rscale, Ystagec));
1281:   PetscCall(MatRestrict(restrct, Zdot, Zdotc));
1282:   PetscCall(VecPointwiseMult(Zdotc, rscale, Zdotc));
1283:   PetscCall(MatRestrict(restrct, Zstage, Zstagec));
1284:   PetscCall(VecPointwiseMult(Zstagec, rscale, Zstagec));
1285:   PetscCall(TSRosWRestoreVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
1286:   PetscCall(TSRosWRestoreVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
1287:   PetscFunctionReturn(PETSC_SUCCESS);
1288: }

1290: static PetscErrorCode DMSubDomainHook_TSRosW(DM fine, DM coarse, void *ctx)
1291: {
1292:   PetscFunctionBegin;
1293:   PetscFunctionReturn(PETSC_SUCCESS);
1294: }

1296: static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
1297: {
1298:   TS  ts = (TS)ctx;
1299:   Vec Ydot, Zdot, Ystage, Zstage;
1300:   Vec Ydots, Zdots, Ystages, Zstages;

1302:   PetscFunctionBegin;
1303:   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
1304:   PetscCall(TSRosWGetVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));

1306:   PetscCall(VecScatterBegin(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
1307:   PetscCall(VecScatterEnd(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));

1309:   PetscCall(VecScatterBegin(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
1310:   PetscCall(VecScatterEnd(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));

1312:   PetscCall(VecScatterBegin(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
1313:   PetscCall(VecScatterEnd(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));

1315:   PetscCall(VecScatterBegin(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
1316:   PetscCall(VecScatterEnd(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));

1318:   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
1319:   PetscCall(TSRosWRestoreVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
1320:   PetscFunctionReturn(PETSC_SUCCESS);
1321: }

1323: /*
1324:   This defines the nonlinear equation that is to be solved with SNES
1325:   G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1326: */
1327: static PetscErrorCode SNESTSFormFunction_RosW(SNES snes, Vec U, Vec F, TS ts)
1328: {
1329:   TS_RosW  *ros = (TS_RosW *)ts->data;
1330:   Vec       Ydot, Zdot, Ystage, Zstage;
1331:   PetscReal shift = ros->scoeff / ts->time_step;
1332:   DM        dm, dmsave;

1334:   PetscFunctionBegin;
1335:   PetscCall(SNESGetDM(snes, &dm));
1336:   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1337:   PetscCall(VecWAXPY(Ydot, shift, U, Zdot));   /* Ydot = shift*U + Zdot */
1338:   PetscCall(VecWAXPY(Ystage, 1.0, U, Zstage)); /* Ystage = U + Zstage */
1339:   dmsave = ts->dm;
1340:   ts->dm = dm;
1341:   PetscCall(TSComputeIFunction(ts, ros->stage_time, Ystage, Ydot, F, PETSC_FALSE));
1342:   ts->dm = dmsave;
1343:   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1344:   PetscFunctionReturn(PETSC_SUCCESS);
1345: }

1347: static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes, Vec U, Mat A, Mat B, TS ts)
1348: {
1349:   TS_RosW  *ros = (TS_RosW *)ts->data;
1350:   Vec       Ydot, Zdot, Ystage, Zstage;
1351:   PetscReal shift = ros->scoeff / ts->time_step;
1352:   DM        dm, dmsave;

1354:   PetscFunctionBegin;
1355:   /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
1356:   PetscCall(SNESGetDM(snes, &dm));
1357:   PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1358:   dmsave = ts->dm;
1359:   ts->dm = dm;
1360:   PetscCall(TSComputeIJacobian(ts, ros->stage_time, Ystage, Ydot, shift, A, B, PETSC_TRUE));
1361:   ts->dm = dmsave;
1362:   PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1363:   PetscFunctionReturn(PETSC_SUCCESS);
1364: }

1366: static PetscErrorCode TSRosWTableauSetUp(TS ts)
1367: {
1368:   TS_RosW    *ros = (TS_RosW *)ts->data;
1369:   RosWTableau tab = ros->tableau;

1371:   PetscFunctionBegin;
1372:   PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ros->Y));
1373:   PetscCall(PetscMalloc1(tab->s, &ros->work));
1374:   PetscFunctionReturn(PETSC_SUCCESS);
1375: }

1377: static PetscErrorCode TSSetUp_RosW(TS ts)
1378: {
1379:   TS_RosW      *ros = (TS_RosW *)ts->data;
1380:   DM            dm;
1381:   SNES          snes;
1382:   TSRHSJacobian rhsjacobian;

1384:   PetscFunctionBegin;
1385:   PetscCall(TSRosWTableauSetUp(ts));
1386:   PetscCall(VecDuplicate(ts->vec_sol, &ros->Ydot));
1387:   PetscCall(VecDuplicate(ts->vec_sol, &ros->Ystage));
1388:   PetscCall(VecDuplicate(ts->vec_sol, &ros->Zdot));
1389:   PetscCall(VecDuplicate(ts->vec_sol, &ros->Zstage));
1390:   PetscCall(VecDuplicate(ts->vec_sol, &ros->vec_sol_prev));
1391:   PetscCall(TSGetDM(ts, &dm));
1392:   PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
1393:   PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1394:   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1395:   PetscCall(TSGetSNES(ts, &snes));
1396:   if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
1397:   PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL));
1398:   if (rhsjacobian == TSComputeRHSJacobianConstant) {
1399:     Mat Amat, Pmat;

1401:     /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */
1402:     PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL));
1403:     if (Amat && Amat == ts->Arhs) {
1404:       if (Amat == Pmat) {
1405:         PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
1406:         PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL));
1407:       } else {
1408:         PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
1409:         PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL));
1410:         if (Pmat && Pmat == ts->Brhs) {
1411:           PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat));
1412:           PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL));
1413:           PetscCall(MatDestroy(&Pmat));
1414:         }
1415:       }
1416:       PetscCall(MatDestroy(&Amat));
1417:     }
1418:   }
1419:   PetscFunctionReturn(PETSC_SUCCESS);
1420: }
1421: /*------------------------------------------------------------*/

1423: static PetscErrorCode TSSetFromOptions_RosW(TS ts, PetscOptionItems *PetscOptionsObject)
1424: {
1425:   TS_RosW *ros = (TS_RosW *)ts->data;
1426:   SNES     snes;

1428:   PetscFunctionBegin;
1429:   PetscOptionsHeadBegin(PetscOptionsObject, "RosW ODE solver options");
1430:   {
1431:     RosWTableauLink link;
1432:     PetscInt        count, choice;
1433:     PetscBool       flg;
1434:     const char    **namelist;

1436:     for (link = RosWTableauList, count = 0; link; link = link->next, count++)
1437:       ;
1438:     PetscCall(PetscMalloc1(count, (char ***)&namelist));
1439:     for (link = RosWTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
1440:     PetscCall(PetscOptionsEList("-ts_rosw_type", "Family of Rosenbrock-W method", "TSRosWSetType", (const char *const *)namelist, count, ros->tableau->name, &choice, &flg));
1441:     if (flg) PetscCall(TSRosWSetType(ts, namelist[choice]));
1442:     PetscCall(PetscFree(namelist));

1444:     PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian", "Recompute the Jacobian at each stage", "TSRosWSetRecomputeJacobian", ros->recompute_jacobian, &ros->recompute_jacobian, NULL));
1445:   }
1446:   PetscOptionsHeadEnd();
1447:   /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1448:   PetscCall(TSGetSNES(ts, &snes));
1449:   if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
1450:   PetscFunctionReturn(PETSC_SUCCESS);
1451: }

1453: static PetscErrorCode TSView_RosW(TS ts, PetscViewer viewer)
1454: {
1455:   TS_RosW  *ros = (TS_RosW *)ts->data;
1456:   PetscBool iascii;

1458:   PetscFunctionBegin;
1459:   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1460:   if (iascii) {
1461:     RosWTableau tab = ros->tableau;
1462:     TSRosWType  rostype;
1463:     char        buf[512];
1464:     PetscInt    i;
1465:     PetscReal   abscissa[512];
1466:     PetscCall(TSRosWGetType(ts, &rostype));
1467:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Rosenbrock-W %s\n", rostype));
1468:     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ASum));
1469:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa of A       = %s\n", buf));
1470:     for (i = 0; i < tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1471:     PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, abscissa));
1472:     PetscCall(PetscViewerASCIIPrintf(viewer, "  Abscissa of A+Gamma = %s\n", buf));
1473:   }
1474:   PetscFunctionReturn(PETSC_SUCCESS);
1475: }

1477: static PetscErrorCode TSLoad_RosW(TS ts, PetscViewer viewer)
1478: {
1479:   SNES    snes;
1480:   TSAdapt adapt;

1482:   PetscFunctionBegin;
1483:   PetscCall(TSGetAdapt(ts, &adapt));
1484:   PetscCall(TSAdaptLoad(adapt, viewer));
1485:   PetscCall(TSGetSNES(ts, &snes));
1486:   PetscCall(SNESLoad(snes, viewer));
1487:   /* function and Jacobian context for SNES when used with TS is always ts object */
1488:   PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
1489:   PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
1490:   PetscFunctionReturn(PETSC_SUCCESS);
1491: }

1493: /*@C
1494:   TSRosWSetType - Set the type of Rosenbrock-W, `TSROSW`, scheme

1496:   Logically Collective

1498:   Input Parameters:
1499: + ts       - timestepping context
1500: - roswtype - type of Rosenbrock-W scheme

1502:   Level: beginner

1504: .seealso: [](ch_ts), `TSRosWGetType()`, `TSROSW`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWARK3`
1505: @*/
1506: PetscErrorCode TSRosWSetType(TS ts, TSRosWType roswtype)
1507: {
1508:   PetscFunctionBegin;
1510:   PetscAssertPointer(roswtype, 2);
1511:   PetscTryMethod(ts, "TSRosWSetType_C", (TS, TSRosWType), (ts, roswtype));
1512:   PetscFunctionReturn(PETSC_SUCCESS);
1513: }

1515: /*@C
1516:   TSRosWGetType - Get the type of Rosenbrock-W scheme

1518:   Logically Collective

1520:   Input Parameter:
1521: . ts - timestepping context

1523:   Output Parameter:
1524: . rostype - type of Rosenbrock-W scheme

1526:   Level: intermediate

1528: .seealso: [](ch_ts), `TSRosWType`, `TSRosWSetType()`
1529: @*/
1530: PetscErrorCode TSRosWGetType(TS ts, TSRosWType *rostype)
1531: {
1532:   PetscFunctionBegin;
1534:   PetscUseMethod(ts, "TSRosWGetType_C", (TS, TSRosWType *), (ts, rostype));
1535:   PetscFunctionReturn(PETSC_SUCCESS);
1536: }

1538: /*@C
1539:   TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.

1541:   Logically Collective

1543:   Input Parameters:
1544: + ts  - timestepping context
1545: - flg - `PETSC_TRUE` to recompute the Jacobian at each stage

1547:   Level: intermediate

1549: .seealso: [](ch_ts), `TSRosWType`, `TSRosWGetType()`
1550: @*/
1551: PetscErrorCode TSRosWSetRecomputeJacobian(TS ts, PetscBool flg)
1552: {
1553:   PetscFunctionBegin;
1555:   PetscTryMethod(ts, "TSRosWSetRecomputeJacobian_C", (TS, PetscBool), (ts, flg));
1556:   PetscFunctionReturn(PETSC_SUCCESS);
1557: }

1559: static PetscErrorCode TSRosWGetType_RosW(TS ts, TSRosWType *rostype)
1560: {
1561:   TS_RosW *ros = (TS_RosW *)ts->data;

1563:   PetscFunctionBegin;
1564:   *rostype = ros->tableau->name;
1565:   PetscFunctionReturn(PETSC_SUCCESS);
1566: }

1568: static PetscErrorCode TSRosWSetType_RosW(TS ts, TSRosWType rostype)
1569: {
1570:   TS_RosW        *ros = (TS_RosW *)ts->data;
1571:   PetscBool       match;
1572:   RosWTableauLink link;

1574:   PetscFunctionBegin;
1575:   if (ros->tableau) {
1576:     PetscCall(PetscStrcmp(ros->tableau->name, rostype, &match));
1577:     if (match) PetscFunctionReturn(PETSC_SUCCESS);
1578:   }
1579:   for (link = RosWTableauList; link; link = link->next) {
1580:     PetscCall(PetscStrcmp(link->tab.name, rostype, &match));
1581:     if (match) {
1582:       if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts));
1583:       ros->tableau = &link->tab;
1584:       if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts));
1585:       ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1586:       PetscFunctionReturn(PETSC_SUCCESS);
1587:     }
1588:   }
1589:   SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rostype);
1590: }

1592: static PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts, PetscBool flg)
1593: {
1594:   TS_RosW *ros = (TS_RosW *)ts->data;

1596:   PetscFunctionBegin;
1597:   ros->recompute_jacobian = flg;
1598:   PetscFunctionReturn(PETSC_SUCCESS);
1599: }

1601: static PetscErrorCode TSDestroy_RosW(TS ts)
1602: {
1603:   PetscFunctionBegin;
1604:   PetscCall(TSReset_RosW(ts));
1605:   if (ts->dm) {
1606:     PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
1607:     PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1608:   }
1609:   PetscCall(PetscFree(ts->data));
1610:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", NULL));
1611:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", NULL));
1612:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", NULL));
1613:   PetscFunctionReturn(PETSC_SUCCESS);
1614: }

1616: /* ------------------------------------------------------------ */
1617: /*MC
1618:       TSROSW - ODE solver using Rosenbrock-W schemes

1620:   These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1621:   nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1622:   of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`.

1624:   Level: beginner

1626:   Notes:
1627:   This method currently only works with autonomous ODE and DAE.

1629:   Consider trying `TSARKIMEX` if the stiff part is strongly nonlinear.

1631:   Since this uses a single linear solve per time-step if you wish to lag the jacobian or preconditioner computation you must use also -snes_lag_jacobian_persists true or -snes_lag_jacobian_preconditioner true

1633:   Developer Notes:
1634:   Rosenbrock-W methods are typically specified for autonomous ODE

1636: $  udot = f(u)

1638:   by the stage equations

1640: $  k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j

1642:   and step completion formula

1644: $  u_1 = u_0 + sum_j b_j k_j

1646:   with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u)
1647:   and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
1648:   we define new variables for the stage equations

1650: $  y_i = gamma_ij k_j

1652:   The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define

1654: $  A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1}

1656:   to rewrite the method as

1658: .vb
1659:   [M/(h gamma_ii) - J] y_i = f(u_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
1660:   u_1 = u_0 + sum_j bt_j y_j
1661: .ve

1663:    where we have introduced the mass matrix M. Continue by defining

1665: $  ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j

1667:    or, more compactly in tensor notation

1669: $  Ydot = 1/h (Gamma^{-1} \otimes I) Y .

1671:    Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
1672:    stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
1673:    equation

1675: $  g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0

1677:    with initial guess y_i = 0.

1679: .seealso: [](ch_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSRosWSetType()`, `TSRosWRegister()`, `TSROSWTHETA1`, `TSROSWTHETA2`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`,
1680:           `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`, `TSType`
1681: M*/
1682: PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1683: {
1684:   TS_RosW *ros;

1686:   PetscFunctionBegin;
1687:   PetscCall(TSRosWInitializePackage());

1689:   ts->ops->reset          = TSReset_RosW;
1690:   ts->ops->destroy        = TSDestroy_RosW;
1691:   ts->ops->view           = TSView_RosW;
1692:   ts->ops->load           = TSLoad_RosW;
1693:   ts->ops->setup          = TSSetUp_RosW;
1694:   ts->ops->step           = TSStep_RosW;
1695:   ts->ops->interpolate    = TSInterpolate_RosW;
1696:   ts->ops->evaluatestep   = TSEvaluateStep_RosW;
1697:   ts->ops->rollback       = TSRollBack_RosW;
1698:   ts->ops->setfromoptions = TSSetFromOptions_RosW;
1699:   ts->ops->snesfunction   = SNESTSFormFunction_RosW;
1700:   ts->ops->snesjacobian   = SNESTSFormJacobian_RosW;

1702:   ts->usessnes = PETSC_TRUE;

1704:   PetscCall(PetscNew(&ros));
1705:   ts->data = (void *)ros;

1707:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", TSRosWGetType_RosW));
1708:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", TSRosWSetType_RosW));
1709:   PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", TSRosWSetRecomputeJacobian_RosW));

1711:   PetscCall(TSRosWSetType(ts, TSRosWDefault));
1712:   PetscFunctionReturn(PETSC_SUCCESS);
1713: }