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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: [](chapter_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: {
330: {0, 0},
331: {1., 0}
332: },
333: Gamma[2][2] = {{1.707106781186547524401, 0}, {-2. * 1.707106781186547524401, 1.707106781186547524401}}, b[2] = {0.5, 0.5}, b1[2] = {1.0, 0.0};
334: PetscReal binterpt[2][2];
335: binterpt[0][0] = 1.707106781186547524401 - 1.0;
336: binterpt[1][0] = 2.0 - 1.707106781186547524401;
337: binterpt[0][1] = 1.707106781186547524401 - 1.5;
338: binterpt[1][1] = 1.5 - 1.707106781186547524401;
340: PetscCall(TSRosWRegister(TSROSW2P, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
341: }
342: {
343: /*const PetscReal g = 1. - 1./PetscSqrtReal(2.0); Direct evaluation: 0.2928932188134524755992. Used for setting up arrays of values known at compile time below. */
344: const PetscReal A[2][2] =
345: {
346: {0, 0},
347: {1., 0}
348: },
349: Gamma[2][2] = {{0.2928932188134524755992, 0}, {-2. * 0.2928932188134524755992, 0.2928932188134524755992}}, b[2] = {0.5, 0.5}, b1[2] = {1.0, 0.0};
350: PetscReal binterpt[2][2];
351: binterpt[0][0] = 0.2928932188134524755992 - 1.0;
352: binterpt[1][0] = 2.0 - 0.2928932188134524755992;
353: binterpt[0][1] = 0.2928932188134524755992 - 1.5;
354: binterpt[1][1] = 1.5 - 0.2928932188134524755992;
356: PetscCall(TSRosWRegister(TSROSW2M, 2, 2, &A[0][0], &Gamma[0][0], b, b1, 2, &binterpt[0][0]));
357: }
358: {
359: /*const PetscReal g = 7.8867513459481287e-01; Directly written in-place below */
360: PetscReal binterpt[3][2];
361: const PetscReal A[3][3] =
362: {
363: {0, 0, 0},
364: {1.5773502691896257e+00, 0, 0},
365: {0.5, 0, 0}
366: },
367: Gamma[3][3] = {{7.8867513459481287e-01, 0, 0}, {-1.5773502691896257e+00, 7.8867513459481287e-01, 0}, {-6.7075317547305480e-01, -1.7075317547305482e-01, 7.8867513459481287e-01}}, b[3] = {1.0566243270259355e-01, 4.9038105676657971e-02, 8.4529946162074843e-01}, b2[3] = {-1.7863279495408180e-01, 1. / 3., 8.4529946162074843e-01};
369: binterpt[0][0] = -0.8094010767585034;
370: binterpt[1][0] = -0.5;
371: binterpt[2][0] = 2.3094010767585034;
372: binterpt[0][1] = 0.9641016151377548;
373: binterpt[1][1] = 0.5;
374: binterpt[2][1] = -1.4641016151377548;
376: PetscCall(TSRosWRegister(TSROSWRA3PW, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
377: }
378: {
379: PetscReal binterpt[4][3];
380: /*const PetscReal g = 4.3586652150845900e-01; Directly written in-place below */
381: const PetscReal A[4][4] =
382: {
383: {0, 0, 0, 0},
384: {8.7173304301691801e-01, 0, 0, 0},
385: {8.4457060015369423e-01, -1.1299064236484185e-01, 0, 0},
386: {0, 0, 1., 0}
387: },
388: Gamma[4][4] = {{4.3586652150845900e-01, 0, 0, 0}, {-8.7173304301691801e-01, 4.3586652150845900e-01, 0, 0}, {-9.0338057013044082e-01, 5.4180672388095326e-02, 4.3586652150845900e-01, 0}, {2.4212380706095346e-01, -1.2232505839045147e+00, 5.4526025533510214e-01, 4.3586652150845900e-01}}, b[4] = {2.4212380706095346e-01, -1.2232505839045147e+00, 1.5452602553351020e+00, 4.3586652150845900e-01}, b2[4] = {3.7810903145819369e-01, -9.6042292212423178e-02, 5.0000000000000000e-01, 2.1793326075422950e-01};
390: binterpt[0][0] = 1.0564298455794094;
391: binterpt[1][0] = 2.296429974281067;
392: binterpt[2][0] = -1.307599564525376;
393: binterpt[3][0] = -1.045260255335102;
394: binterpt[0][1] = -1.3864882699759573;
395: binterpt[1][1] = -8.262611700275677;
396: binterpt[2][1] = 7.250979895056055;
397: binterpt[3][1] = 2.398120075195581;
398: binterpt[0][2] = 0.5721822314575016;
399: binterpt[1][2] = 4.742931142090097;
400: binterpt[2][2] = -4.398120075195578;
401: binterpt[3][2] = -0.9169932983520199;
403: PetscCall(TSRosWRegister(TSROSWRA34PW2, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
404: }
405: {
406: /* const PetscReal g = 0.5; Directly written in-place below */
407: const PetscReal A[4][4] =
408: {
409: {0, 0, 0, 0},
410: {0, 0, 0, 0},
411: {1., 0, 0, 0},
412: {0.75, -0.25, 0.5, 0}
413: },
414: Gamma[4][4] = {{0.5, 0, 0, 0}, {1., 0.5, 0, 0}, {-0.25, -0.25, 0.5, 0}, {1. / 12, 1. / 12, -2. / 3, 0.5}}, b[4] = {5. / 6, -1. / 6, -1. / 6, 0.5}, b2[4] = {0.75, -0.25, 0.5, 0};
416: PetscCall(TSRosWRegister(TSROSWRODAS3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 0, NULL));
417: }
418: {
419: /*const PetscReal g = 0.43586652150845899941601945119356; Directly written in-place below */
420: const PetscReal A[3][3] =
421: {
422: {0, 0, 0},
423: {0.43586652150845899941601945119356, 0, 0},
424: {0.43586652150845899941601945119356, 0, 0}
425: },
426: Gamma[3][3] = {{0.43586652150845899941601945119356, 0, 0}, {-0.19294655696029095575009695436041, 0.43586652150845899941601945119356, 0}, {0, 1.74927148125794685173529749738960, 0.43586652150845899941601945119356}}, b[3] = {-0.75457412385404315829818998646589, 1.94100407061964420292840123379419, -0.18642994676560104463021124732829}, b2[3] = {-1.53358745784149585370766523913002, 2.81745131148625772213931745457622, -0.28386385364476186843165221544619};
428: PetscReal binterpt[3][2];
429: binterpt[0][0] = 3.793692883777660870425141387941;
430: binterpt[1][0] = -2.918692883777660870425141387941;
431: binterpt[2][0] = 0.125;
432: binterpt[0][1] = -0.725741064379812106687651020584;
433: binterpt[1][1] = 0.559074397713145440020984353917;
434: binterpt[2][1] = 0.16666666666666666666666666666667;
436: PetscCall(TSRosWRegister(TSROSWSANDU3, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
437: }
438: {
439: /*const PetscReal s3 = PetscSqrtReal(3.),g = (3.0+s3)/6.0;
440: * Direct evaluation: s3 = 1.732050807568877293527;
441: * g = 0.7886751345948128822546;
442: * Values are directly inserted below to ensure availability at compile time (compiler warnings otherwise...) */
443: const PetscReal A[3][3] =
444: {
445: {0, 0, 0},
446: {1, 0, 0},
447: {0.25, 0.25, 0}
448: },
449: Gamma[3][3] = {{0, 0, 0}, {(-3.0 - 1.732050807568877293527) / 6.0, 0.7886751345948128822546, 0}, {(-3.0 - 1.732050807568877293527) / 24.0, (-3.0 - 1.732050807568877293527) / 8.0, 0.7886751345948128822546}}, b[3] = {1. / 6., 1. / 6., 2. / 3.}, b2[3] = {1. / 4., 1. / 4., 1. / 2.};
450: PetscReal binterpt[3][2];
452: binterpt[0][0] = 0.089316397477040902157517886164709;
453: binterpt[1][0] = -0.91068360252295909784248211383529;
454: binterpt[2][0] = 1.8213672050459181956849642276706;
455: binterpt[0][1] = 0.077350269189625764509148780501957;
456: binterpt[1][1] = 1.077350269189625764509148780502;
457: binterpt[2][1] = -1.1547005383792515290182975610039;
459: PetscCall(TSRosWRegister(TSROSWASSP3P3S1C, 3, 3, &A[0][0], &Gamma[0][0], b, b2, 2, &binterpt[0][0]));
460: }
462: {
463: const PetscReal A[4][4] =
464: {
465: {0, 0, 0, 0},
466: {1. / 2., 0, 0, 0},
467: {1. / 2., 1. / 2., 0, 0},
468: {1. / 6., 1. / 6., 1. / 6., 0}
469: },
470: Gamma[4][4] = {{1. / 2., 0, 0, 0}, {0.0, 1. / 4., 0, 0}, {-2., -2. / 3., 2. / 3., 0}, {1. / 2., 5. / 36., -2. / 9, 0}}, b[4] = {1. / 6., 1. / 6., 1. / 6., 1. / 2.}, b2[4] = {1. / 8., 3. / 4., 1. / 8., 0};
471: PetscReal binterpt[4][3];
473: binterpt[0][0] = 6.25;
474: binterpt[1][0] = -30.25;
475: binterpt[2][0] = 1.75;
476: binterpt[3][0] = 23.25;
477: binterpt[0][1] = -9.75;
478: binterpt[1][1] = 58.75;
479: binterpt[2][1] = -3.25;
480: binterpt[3][1] = -45.75;
481: binterpt[0][2] = 3.6666666666666666666666666666667;
482: binterpt[1][2] = -28.333333333333333333333333333333;
483: binterpt[2][2] = 1.6666666666666666666666666666667;
484: binterpt[3][2] = 23.;
486: PetscCall(TSRosWRegister(TSROSWLASSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
487: }
489: {
490: const PetscReal A[4][4] =
491: {
492: {0, 0, 0, 0},
493: {1. / 2., 0, 0, 0},
494: {1. / 2., 1. / 2., 0, 0},
495: {1. / 6., 1. / 6., 1. / 6., 0}
496: },
497: Gamma[4][4] = {{1. / 2., 0, 0, 0}, {0.0, 3. / 4., 0, 0}, {-2. / 3., -23. / 9., 2. / 9., 0}, {1. / 18., 65. / 108., -2. / 27, 0}}, b[4] = {1. / 6., 1. / 6., 1. / 6., 1. / 2.}, b2[4] = {3. / 16., 10. / 16., 3. / 16., 0};
498: PetscReal binterpt[4][3];
500: binterpt[0][0] = 1.6911764705882352941176470588235;
501: binterpt[1][0] = 3.6813725490196078431372549019608;
502: binterpt[2][0] = 0.23039215686274509803921568627451;
503: binterpt[3][0] = -4.6029411764705882352941176470588;
504: binterpt[0][1] = -0.95588235294117647058823529411765;
505: binterpt[1][1] = -6.2401960784313725490196078431373;
506: binterpt[2][1] = -0.31862745098039215686274509803922;
507: binterpt[3][1] = 7.5147058823529411764705882352941;
508: binterpt[0][2] = -0.56862745098039215686274509803922;
509: binterpt[1][2] = 2.7254901960784313725490196078431;
510: binterpt[2][2] = 0.25490196078431372549019607843137;
511: binterpt[3][2] = -2.4117647058823529411764705882353;
513: PetscCall(TSRosWRegister(TSROSWLLSSP3P4S2C, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
514: }
516: {
517: PetscReal A[4][4], Gamma[4][4], b[4], b2[4];
518: PetscReal binterpt[4][3];
520: Gamma[0][0] = 0.4358665215084589994160194475295062513822671686978816;
521: Gamma[0][1] = 0;
522: Gamma[0][2] = 0;
523: Gamma[0][3] = 0;
524: Gamma[1][0] = -1.997527830934941248426324674704153457289527280554476;
525: Gamma[1][1] = 0.4358665215084589994160194475295062513822671686978816;
526: Gamma[1][2] = 0;
527: Gamma[1][3] = 0;
528: Gamma[2][0] = -1.007948511795029620852002345345404191008352770119903;
529: Gamma[2][1] = -0.004648958462629345562774289390054679806993396798458131;
530: Gamma[2][2] = 0.4358665215084589994160194475295062513822671686978816;
531: Gamma[2][3] = 0;
532: Gamma[3][0] = -0.6685429734233467180451604600279552604364311322650783;
533: Gamma[3][1] = 0.6056625986449338476089525334450053439525178740492984;
534: Gamma[3][2] = -0.9717899277217721234705114616271378792182450260943198;
535: Gamma[3][3] = 0;
537: A[0][0] = 0;
538: A[0][1] = 0;
539: A[0][2] = 0;
540: A[0][3] = 0;
541: A[1][0] = 0.8717330430169179988320388950590125027645343373957631;
542: A[1][1] = 0;
543: A[1][2] = 0;
544: A[1][3] = 0;
545: A[2][0] = 0.5275890119763004115618079766722914408876108660811028;
546: A[2][1] = 0.07241098802369958843819203208518599088698057726988732;
547: A[2][2] = 0;
548: A[2][3] = 0;
549: A[3][0] = 0.3990960076760701320627260685975778145384666450351314;
550: A[3][1] = -0.4375576546135194437228463747348862825846903771419953;
551: A[3][2] = 1.038461646937449311660120300601880176655352737312713;
552: A[3][3] = 0;
554: b[0] = 0.1876410243467238251612921333138006734899663569186926;
555: b[1] = -0.5952974735769549480478230473706443582188442040780541;
556: b[2] = 0.9717899277217721234705114616271378792182450260943198;
557: b[3] = 0.4358665215084589994160194475295062513822671686978816;
559: b2[0] = 0.2147402862233891404862383521089097657790734483804460;
560: b2[1] = -0.4851622638849390928209050538171743017757490232519684;
561: b2[2] = 0.8687250025203875511662123688667549217531982787600080;
562: b2[3] = 0.4016969751411624011684543450940068201770721128357014;
564: binterpt[0][0] = 2.2565812720167954547104627844105;
565: binterpt[1][0] = 1.349166413351089573796243820819;
566: binterpt[2][0] = -2.4695174540533503758652847586647;
567: binterpt[3][0] = -0.13623023131453465264142184656474;
568: binterpt[0][1] = -3.0826699111559187902922463354557;
569: binterpt[1][1] = -2.4689115685996042534544925650515;
570: binterpt[2][1] = 5.7428279814696677152129332773553;
571: binterpt[3][1] = -0.19124650171414467146619437684812;
572: binterpt[0][2] = 1.0137296634858471607430756831148;
573: binterpt[1][2] = 0.52444768167155973161042570784064;
574: binterpt[2][2] = -2.3015205996945452158771370439586;
575: binterpt[3][2] = 0.76334325453713832352363565300308;
577: PetscCall(TSRosWRegister(TSROSWARK3, 3, 4, &A[0][0], &Gamma[0][0], b, b2, 3, &binterpt[0][0]));
578: }
579: PetscCall(TSRosWRegisterRos4(TSROSWGRK4T, 0.231, PETSC_DEFAULT, PETSC_DEFAULT, 0, -0.1282612945269037e+01));
580: PetscCall(TSRosWRegisterRos4(TSROSWSHAMP4, 0.5, PETSC_DEFAULT, PETSC_DEFAULT, 0, 125. / 108.));
581: PetscCall(TSRosWRegisterRos4(TSROSWVELDD4, 0.22570811482256823492, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.355958941201148));
582: PetscCall(TSRosWRegisterRos4(TSROSW4L, 0.57282, PETSC_DEFAULT, PETSC_DEFAULT, 0, -1.093502252409163));
583: PetscFunctionReturn(PETSC_SUCCESS);
584: }
586: /*@C
587: TSRosWRegisterDestroy - Frees the list of schemes that were registered by `TSRosWRegister()`.
589: Not Collective
591: Level: advanced
593: .seealso: [](chapter_ts), `TSRosWRegister()`, `TSRosWRegisterAll()`
594: @*/
595: PetscErrorCode TSRosWRegisterDestroy(void)
596: {
597: RosWTableauLink link;
599: PetscFunctionBegin;
600: while ((link = RosWTableauList)) {
601: RosWTableau t = &link->tab;
602: RosWTableauList = link->next;
603: PetscCall(PetscFree5(t->A, t->Gamma, t->b, t->ASum, t->GammaSum));
604: PetscCall(PetscFree5(t->At, t->bt, t->GammaInv, t->GammaZeroDiag, t->GammaExplicitCorr));
605: PetscCall(PetscFree2(t->bembed, t->bembedt));
606: PetscCall(PetscFree(t->binterpt));
607: PetscCall(PetscFree(t->name));
608: PetscCall(PetscFree(link));
609: }
610: TSRosWRegisterAllCalled = PETSC_FALSE;
611: PetscFunctionReturn(PETSC_SUCCESS);
612: }
614: /*@C
615: TSRosWInitializePackage - This function initializes everything in the `TSROSW` package. It is called
616: from `TSInitializePackage()`.
618: Level: developer
620: .seealso: [](chapter_ts), `TSROSW`, `PetscInitialize()`, `TSRosWFinalizePackage()`
621: @*/
622: PetscErrorCode TSRosWInitializePackage(void)
623: {
624: PetscFunctionBegin;
625: if (TSRosWPackageInitialized) PetscFunctionReturn(PETSC_SUCCESS);
626: TSRosWPackageInitialized = PETSC_TRUE;
627: PetscCall(TSRosWRegisterAll());
628: PetscCall(PetscRegisterFinalize(TSRosWFinalizePackage));
629: PetscFunctionReturn(PETSC_SUCCESS);
630: }
632: /*@C
633: TSRosWFinalizePackage - This function destroys everything in the `TSROSW` package. It is
634: called from `PetscFinalize()`.
636: Level: developer
638: .seealso: [](chapter_ts), `TSROSW`, `PetscFinalize()`, `TSRosWInitializePackage()`
639: @*/
640: PetscErrorCode TSRosWFinalizePackage(void)
641: {
642: PetscFunctionBegin;
643: TSRosWPackageInitialized = PETSC_FALSE;
644: PetscCall(TSRosWRegisterDestroy());
645: PetscFunctionReturn(PETSC_SUCCESS);
646: }
648: /*@C
649: TSRosWRegister - register a `TSROSW`, Rosenbrock W scheme by providing the entries in the Butcher tableau and optionally embedded approximations and interpolation
651: Not Collective, but the same schemes should be registered on all processes on which they will be used
653: Input Parameters:
654: + name - identifier for method
655: . order - approximation order of method
656: . s - number of stages, this is the dimension of the matrices below
657: . A - Table of propagated stage coefficients (dimension s*s, row-major), strictly lower triangular
658: . Gamma - Table of coefficients in implicit stage equations (dimension s*s, row-major), lower triangular with nonzero diagonal
659: . b - Step completion table (dimension s)
660: . bembed - Step completion table for a scheme of order one less (dimension s, NULL if no embedded scheme is available)
661: . pinterp - Order of the interpolation scheme, equal to the number of columns of binterpt
662: - binterpt - Coefficients of the interpolation formula (dimension s*pinterp)
664: Level: advanced
666: Note:
667: Several Rosenbrock W methods are provided, this function is only needed to create new methods.
669: .seealso: [](chapter_ts), `TSROSW`
670: @*/
671: 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[])
672: {
673: RosWTableauLink link;
674: RosWTableau t;
675: PetscInt i, j, k;
676: PetscScalar *GammaInv;
678: PetscFunctionBegin;
685: PetscCall(TSRosWInitializePackage());
686: PetscCall(PetscNew(&link));
687: t = &link->tab;
688: PetscCall(PetscStrallocpy(name, &t->name));
689: t->order = order;
690: t->s = s;
691: PetscCall(PetscMalloc5(s * s, &t->A, s * s, &t->Gamma, s, &t->b, s, &t->ASum, s, &t->GammaSum));
692: PetscCall(PetscMalloc5(s * s, &t->At, s, &t->bt, s * s, &t->GammaInv, s, &t->GammaZeroDiag, s * s, &t->GammaExplicitCorr));
693: PetscCall(PetscArraycpy(t->A, A, s * s));
694: PetscCall(PetscArraycpy(t->Gamma, Gamma, s * s));
695: PetscCall(PetscArraycpy(t->GammaExplicitCorr, Gamma, s * s));
696: PetscCall(PetscArraycpy(t->b, b, s));
697: if (bembed) {
698: PetscCall(PetscMalloc2(s, &t->bembed, s, &t->bembedt));
699: PetscCall(PetscArraycpy(t->bembed, bembed, s));
700: }
701: for (i = 0; i < s; i++) {
702: t->ASum[i] = 0;
703: t->GammaSum[i] = 0;
704: for (j = 0; j < s; j++) {
705: t->ASum[i] += A[i * s + j];
706: t->GammaSum[i] += Gamma[i * s + j];
707: }
708: }
709: PetscCall(PetscMalloc1(s * s, &GammaInv)); /* Need to use Scalar for inverse, then convert back to Real */
710: for (i = 0; i < s * s; i++) GammaInv[i] = Gamma[i];
711: for (i = 0; i < s; i++) {
712: if (Gamma[i * s + i] == 0.0) {
713: GammaInv[i * s + i] = 1.0;
714: t->GammaZeroDiag[i] = PETSC_TRUE;
715: } else {
716: t->GammaZeroDiag[i] = PETSC_FALSE;
717: }
718: }
720: switch (s) {
721: case 1:
722: GammaInv[0] = 1. / GammaInv[0];
723: break;
724: case 2:
725: PetscCall(PetscKernel_A_gets_inverse_A_2(GammaInv, 0, PETSC_FALSE, NULL));
726: break;
727: case 3:
728: PetscCall(PetscKernel_A_gets_inverse_A_3(GammaInv, 0, PETSC_FALSE, NULL));
729: break;
730: case 4:
731: PetscCall(PetscKernel_A_gets_inverse_A_4(GammaInv, 0, PETSC_FALSE, NULL));
732: break;
733: case 5: {
734: PetscInt ipvt5[5];
735: MatScalar work5[5 * 5];
736: PetscCall(PetscKernel_A_gets_inverse_A_5(GammaInv, ipvt5, work5, 0, PETSC_FALSE, NULL));
737: break;
738: }
739: case 6:
740: PetscCall(PetscKernel_A_gets_inverse_A_6(GammaInv, 0, PETSC_FALSE, NULL));
741: break;
742: case 7:
743: PetscCall(PetscKernel_A_gets_inverse_A_7(GammaInv, 0, PETSC_FALSE, NULL));
744: break;
745: default:
746: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Not implemented for %" PetscInt_FMT " stages", s);
747: }
748: for (i = 0; i < s * s; i++) t->GammaInv[i] = PetscRealPart(GammaInv[i]);
749: PetscCall(PetscFree(GammaInv));
751: for (i = 0; i < s; i++) {
752: for (k = 0; k < i + 1; k++) {
753: t->GammaExplicitCorr[i * s + k] = (t->GammaExplicitCorr[i * s + k]) * (t->GammaInv[k * s + k]);
754: for (j = k + 1; j < i + 1; j++) t->GammaExplicitCorr[i * s + k] += (t->GammaExplicitCorr[i * s + j]) * (t->GammaInv[j * s + k]);
755: }
756: }
758: for (i = 0; i < s; i++) {
759: for (j = 0; j < s; j++) {
760: t->At[i * s + j] = 0;
761: for (k = 0; k < s; k++) t->At[i * s + j] += t->A[i * s + k] * t->GammaInv[k * s + j];
762: }
763: t->bt[i] = 0;
764: for (j = 0; j < s; j++) t->bt[i] += t->b[j] * t->GammaInv[j * s + i];
765: if (bembed) {
766: t->bembedt[i] = 0;
767: for (j = 0; j < s; j++) t->bembedt[i] += t->bembed[j] * t->GammaInv[j * s + i];
768: }
769: }
770: t->ccfl = 1.0; /* Fix this */
772: t->pinterp = pinterp;
773: PetscCall(PetscMalloc1(s * pinterp, &t->binterpt));
774: PetscCall(PetscArraycpy(t->binterpt, binterpt, s * pinterp));
775: link->next = RosWTableauList;
776: RosWTableauList = link;
777: PetscFunctionReturn(PETSC_SUCCESS);
778: }
780: /*@C
781: TSRosWRegisterRos4 - register a fourth order Rosenbrock scheme by providing parameter choices
783: Not Collective, but the same schemes should be registered on all processes on which they will be used
785: Input Parameters:
786: + name - identifier for method
787: . gamma - leading coefficient (diagonal entry)
788: . a2 - design parameter, see Table 7.2 of Hairer&Wanner
789: . a3 - design parameter or PETSC_DEFAULT to satisfy one of the order five conditions (Eq 7.22)
790: . b3 - design parameter, see Table 7.2 of Hairer&Wanner
791: . beta43 - design parameter or PETSC_DEFAULT to use Equation 7.21 of Hairer&Wanner
792: - e4 - design parameter for embedded method, see coefficient E4 in ros4.f code from Hairer
794: Level: developer
796: Notes:
797: This routine encodes the design of fourth order Rosenbrock methods as described in Hairer and Wanner volume 2.
798: It is used here to implement several methods from the book and can be used to experiment with new methods.
799: It was written this way instead of by copying coefficients in order to provide better than double precision satisfaction of the order conditions.
801: .seealso: [](chapter_ts), `TSRosW`, `TSRosWRegister()`
802: @*/
803: PetscErrorCode TSRosWRegisterRos4(TSRosWType name, PetscReal gamma, PetscReal a2, PetscReal a3, PetscReal b3, PetscReal e4)
804: {
805: /* Declare numeric constants so they can be quad precision without being truncated at double */
806: 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;
807: PetscReal a4, a32, a42, a43, b1, b2, b4, beta2p, beta3p, beta4p, beta32, beta42, beta43, beta32beta2p, beta4jbetajp;
808: PetscReal A[4][4], Gamma[4][4], b[4], bm[4];
809: PetscScalar M[3][3], rhs[3];
811: PetscFunctionBegin;
812: /* Step 1: choose Gamma (input) */
813: /* Step 2: choose a2,a3,a4; b1,b2,b3,b4 to satisfy order conditions */
814: if (a3 == PETSC_DEFAULT) a3 = (one / five - a2 / four) / (one / four - a2 / three); /* Eq 7.22 */
815: a4 = a3; /* consequence of 7.20 */
817: /* Solve order conditions 7.15a, 7.15c, 7.15e */
818: M[0][0] = one;
819: M[0][1] = one;
820: M[0][2] = one; /* 7.15a */
821: M[1][0] = 0.0;
822: M[1][1] = a2 * a2;
823: M[1][2] = a4 * a4; /* 7.15c */
824: M[2][0] = 0.0;
825: M[2][1] = a2 * a2 * a2;
826: M[2][2] = a4 * a4 * a4; /* 7.15e */
827: rhs[0] = one - b3;
828: rhs[1] = one / three - a3 * a3 * b3;
829: rhs[2] = one / four - a3 * a3 * a3 * b3;
830: PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
831: b1 = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
832: b2 = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
833: b4 = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);
835: /* Step 3 */
836: beta43 = (p56 - a2 * p43) / (b4 * a3 * a3 * (a3 - a2)); /* 7.21 */
837: beta32beta2p = p44 / (b4 * beta43); /* 7.15h */
838: beta4jbetajp = (p32 - b3 * beta32beta2p) / b4;
839: M[0][0] = b2;
840: M[0][1] = b3;
841: M[0][2] = b4;
842: M[1][0] = a4 * a4 * beta32beta2p - a3 * a3 * beta4jbetajp;
843: M[1][1] = a2 * a2 * beta4jbetajp;
844: M[1][2] = -a2 * a2 * beta32beta2p;
845: M[2][0] = b4 * beta43 * a3 * a3 - p43;
846: M[2][1] = -b4 * beta43 * a2 * a2;
847: M[2][2] = 0;
848: rhs[0] = one / two - gamma;
849: rhs[1] = 0;
850: rhs[2] = -a2 * a2 * p32;
851: PetscCall(PetscKernel_A_gets_inverse_A_3(&M[0][0], 0, PETSC_FALSE, NULL));
852: beta2p = PetscRealPart(M[0][0] * rhs[0] + M[0][1] * rhs[1] + M[0][2] * rhs[2]);
853: beta3p = PetscRealPart(M[1][0] * rhs[0] + M[1][1] * rhs[1] + M[1][2] * rhs[2]);
854: beta4p = PetscRealPart(M[2][0] * rhs[0] + M[2][1] * rhs[1] + M[2][2] * rhs[2]);
856: /* Step 4: back-substitute */
857: beta32 = beta32beta2p / beta2p;
858: beta42 = (beta4jbetajp - beta43 * beta3p) / beta2p;
860: /* Step 5: 7.15f and 7.20, then 7.16 */
861: a43 = 0;
862: a32 = p42 / (b3 * a3 * beta2p + b4 * a4 * beta2p);
863: a42 = a32;
865: A[0][0] = 0;
866: A[0][1] = 0;
867: A[0][2] = 0;
868: A[0][3] = 0;
869: A[1][0] = a2;
870: A[1][1] = 0;
871: A[1][2] = 0;
872: A[1][3] = 0;
873: A[2][0] = a3 - a32;
874: A[2][1] = a32;
875: A[2][2] = 0;
876: A[2][3] = 0;
877: A[3][0] = a4 - a43 - a42;
878: A[3][1] = a42;
879: A[3][2] = a43;
880: A[3][3] = 0;
881: Gamma[0][0] = gamma;
882: Gamma[0][1] = 0;
883: Gamma[0][2] = 0;
884: Gamma[0][3] = 0;
885: Gamma[1][0] = beta2p - A[1][0];
886: Gamma[1][1] = gamma;
887: Gamma[1][2] = 0;
888: Gamma[1][3] = 0;
889: Gamma[2][0] = beta3p - beta32 - A[2][0];
890: Gamma[2][1] = beta32 - A[2][1];
891: Gamma[2][2] = gamma;
892: Gamma[2][3] = 0;
893: Gamma[3][0] = beta4p - beta42 - beta43 - A[3][0];
894: Gamma[3][1] = beta42 - A[3][1];
895: Gamma[3][2] = beta43 - A[3][2];
896: Gamma[3][3] = gamma;
897: b[0] = b1;
898: b[1] = b2;
899: b[2] = b3;
900: b[3] = b4;
902: /* Construct embedded formula using given e4. We are solving Equation 7.18. */
903: bm[3] = b[3] - e4 * gamma; /* using definition of E4 */
904: bm[2] = (p32 - beta4jbetajp * bm[3]) / (beta32 * beta2p); /* fourth row of 7.18 */
905: bm[1] = (one / two - gamma - beta3p * bm[2] - beta4p * bm[3]) / beta2p; /* second row */
906: bm[0] = one - bm[1] - bm[2] - bm[3]; /* first row */
908: {
909: const PetscReal misfit = a2 * a2 * bm[1] + a3 * a3 * bm[2] + a4 * a4 * bm[3] - one / three;
910: PetscCheck(PetscAbs(misfit) <= PETSC_SMALL, PETSC_COMM_SELF, PETSC_ERR_SUP, "Assumptions violated, could not construct a third order embedded method");
911: }
912: PetscCall(TSRosWRegister(name, 4, 4, &A[0][0], &Gamma[0][0], b, bm, 0, NULL));
913: PetscFunctionReturn(PETSC_SUCCESS);
914: }
916: /*
917: The step completion formula is
919: x1 = x0 + b^T Y
921: where Y is the multi-vector of stages corrections. This function can be called before or after ts->vec_sol has been
922: updated. Suppose we have a completion formula b and an embedded formula be of different order. We can write
924: x1e = x0 + be^T Y
925: = x1 - b^T Y + be^T Y
926: = x1 + (be - b)^T Y
928: so we can evaluate the method of different order even after the step has been optimistically completed.
929: */
930: static PetscErrorCode TSEvaluateStep_RosW(TS ts, PetscInt order, Vec U, PetscBool *done)
931: {
932: TS_RosW *ros = (TS_RosW *)ts->data;
933: RosWTableau tab = ros->tableau;
934: PetscScalar *w = ros->work;
935: PetscInt i;
937: PetscFunctionBegin;
938: if (order == tab->order) {
939: if (ros->status == TS_STEP_INCOMPLETE) { /* Use standard completion formula */
940: PetscCall(VecCopy(ts->vec_sol, U));
941: for (i = 0; i < tab->s; i++) w[i] = tab->bt[i];
942: PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
943: } else PetscCall(VecCopy(ts->vec_sol, U));
944: if (done) *done = PETSC_TRUE;
945: PetscFunctionReturn(PETSC_SUCCESS);
946: } else if (order == tab->order - 1) {
947: if (!tab->bembedt) goto unavailable;
948: if (ros->status == TS_STEP_INCOMPLETE) { /* Use embedded completion formula */
949: PetscCall(VecCopy(ts->vec_sol, U));
950: for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i];
951: PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
952: } else { /* Use rollback-and-recomplete formula (bembedt - bt) */
953: for (i = 0; i < tab->s; i++) w[i] = tab->bembedt[i] - tab->bt[i];
954: PetscCall(VecCopy(ts->vec_sol, U));
955: PetscCall(VecMAXPY(U, tab->s, w, ros->Y));
956: }
957: if (done) *done = PETSC_TRUE;
958: PetscFunctionReturn(PETSC_SUCCESS);
959: }
960: unavailable:
961: if (done) *done = PETSC_FALSE;
962: else
963: 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,
964: tab->order, order);
965: PetscFunctionReturn(PETSC_SUCCESS);
966: }
968: static PetscErrorCode TSRollBack_RosW(TS ts)
969: {
970: TS_RosW *ros = (TS_RosW *)ts->data;
972: PetscFunctionBegin;
973: PetscCall(VecCopy(ros->vec_sol_prev, ts->vec_sol));
974: PetscFunctionReturn(PETSC_SUCCESS);
975: }
977: static PetscErrorCode TSStep_RosW(TS ts)
978: {
979: TS_RosW *ros = (TS_RosW *)ts->data;
980: RosWTableau tab = ros->tableau;
981: const PetscInt s = tab->s;
982: const PetscReal *At = tab->At, *Gamma = tab->Gamma, *ASum = tab->ASum, *GammaInv = tab->GammaInv;
983: const PetscReal *GammaExplicitCorr = tab->GammaExplicitCorr;
984: const PetscBool *GammaZeroDiag = tab->GammaZeroDiag;
985: PetscScalar *w = ros->work;
986: Vec *Y = ros->Y, Ydot = ros->Ydot, Zdot = ros->Zdot, Zstage = ros->Zstage;
987: SNES snes;
988: TSAdapt adapt;
989: PetscInt i, j, its, lits;
990: PetscInt rejections = 0;
991: PetscBool stageok, accept = PETSC_TRUE;
992: PetscReal next_time_step = ts->time_step;
993: PetscInt lag;
995: PetscFunctionBegin;
996: if (!ts->steprollback) PetscCall(VecCopy(ts->vec_sol, ros->vec_sol_prev));
998: ros->status = TS_STEP_INCOMPLETE;
999: while (!ts->reason && ros->status != TS_STEP_COMPLETE) {
1000: const PetscReal h = ts->time_step;
1001: for (i = 0; i < s; i++) {
1002: ros->stage_time = ts->ptime + h * ASum[i];
1003: PetscCall(TSPreStage(ts, ros->stage_time));
1004: if (GammaZeroDiag[i]) {
1005: ros->stage_explicit = PETSC_TRUE;
1006: ros->scoeff = 1.;
1007: } else {
1008: ros->stage_explicit = PETSC_FALSE;
1009: ros->scoeff = 1. / Gamma[i * s + i];
1010: }
1012: PetscCall(VecCopy(ts->vec_sol, Zstage));
1013: for (j = 0; j < i; j++) w[j] = At[i * s + j];
1014: PetscCall(VecMAXPY(Zstage, i, w, Y));
1016: for (j = 0; j < i; j++) w[j] = 1. / h * GammaInv[i * s + j];
1017: PetscCall(VecZeroEntries(Zdot));
1018: PetscCall(VecMAXPY(Zdot, i, w, Y));
1020: /* Initial guess taken from last stage */
1021: PetscCall(VecZeroEntries(Y[i]));
1023: if (!ros->stage_explicit) {
1024: PetscCall(TSGetSNES(ts, &snes));
1025: if (!ros->recompute_jacobian && !i) {
1026: PetscCall(SNESGetLagJacobian(snes, &lag));
1027: if (lag == 1) { /* use did not set a nontrivial lag, so lag over all stages */
1028: PetscCall(SNESSetLagJacobian(snes, -2)); /* Recompute the Jacobian on this solve, but not again for the rest of the stages */
1029: }
1030: }
1031: PetscCall(SNESSolve(snes, NULL, Y[i]));
1032: 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 */ }
1033: PetscCall(SNESGetIterationNumber(snes, &its));
1034: PetscCall(SNESGetLinearSolveIterations(snes, &lits));
1035: ts->snes_its += its;
1036: ts->ksp_its += lits;
1037: } else {
1038: Mat J, Jp;
1039: PetscCall(VecZeroEntries(Ydot)); /* Evaluate Y[i]=G(t,Ydot=0,Zstage) */
1040: PetscCall(TSComputeIFunction(ts, ros->stage_time, Zstage, Ydot, Y[i], PETSC_FALSE));
1041: PetscCall(VecScale(Y[i], -1.0));
1042: PetscCall(VecAXPY(Y[i], -1.0, Zdot)); /*Y[i] = F(Zstage)-Zdot[=GammaInv*Y]*/
1044: PetscCall(VecZeroEntries(Zstage)); /* Zstage = GammaExplicitCorr[i,j] * Y[j] */
1045: for (j = 0; j < i; j++) w[j] = GammaExplicitCorr[i * s + j];
1046: PetscCall(VecMAXPY(Zstage, i, w, Y));
1048: /* Y[i] = Y[i] + Jac*Zstage[=Jac*GammaExplicitCorr[i,j] * Y[j]] */
1049: PetscCall(TSGetIJacobian(ts, &J, &Jp, NULL, NULL));
1050: PetscCall(TSComputeIJacobian(ts, ros->stage_time, ts->vec_sol, Ydot, 0, J, Jp, PETSC_FALSE));
1051: PetscCall(MatMult(J, Zstage, Zdot));
1052: PetscCall(VecAXPY(Y[i], -1.0, Zdot));
1053: ts->ksp_its += 1;
1055: PetscCall(VecScale(Y[i], h));
1056: }
1057: PetscCall(TSPostStage(ts, ros->stage_time, i, Y));
1058: PetscCall(TSGetAdapt(ts, &adapt));
1059: PetscCall(TSAdaptCheckStage(adapt, ts, ros->stage_time, Y[i], &stageok));
1060: if (!stageok) goto reject_step;
1061: }
1063: ros->status = TS_STEP_INCOMPLETE;
1064: PetscCall(TSEvaluateStep_RosW(ts, tab->order, ts->vec_sol, NULL));
1065: ros->status = TS_STEP_PENDING;
1066: PetscCall(TSGetAdapt(ts, &adapt));
1067: PetscCall(TSAdaptCandidatesClear(adapt));
1068: PetscCall(TSAdaptCandidateAdd(adapt, tab->name, tab->order, 1, tab->ccfl, (PetscReal)tab->s, PETSC_TRUE));
1069: PetscCall(TSAdaptChoose(adapt, ts, ts->time_step, NULL, &next_time_step, &accept));
1070: ros->status = accept ? TS_STEP_COMPLETE : TS_STEP_INCOMPLETE;
1071: if (!accept) { /* Roll back the current step */
1072: PetscCall(TSRollBack_RosW(ts));
1073: ts->time_step = next_time_step;
1074: goto reject_step;
1075: }
1077: ts->ptime += ts->time_step;
1078: ts->time_step = next_time_step;
1079: break;
1081: reject_step:
1082: ts->reject++;
1083: accept = PETSC_FALSE;
1084: if (!ts->reason && ++rejections > ts->max_reject && ts->max_reject >= 0) {
1085: ts->reason = TS_DIVERGED_STEP_REJECTED;
1086: PetscCall(PetscInfo(ts, "Step=%" PetscInt_FMT ", step rejections %" PetscInt_FMT " greater than current TS allowed, stopping solve\n", ts->steps, rejections));
1087: }
1088: }
1089: PetscFunctionReturn(PETSC_SUCCESS);
1090: }
1092: static PetscErrorCode TSInterpolate_RosW(TS ts, PetscReal itime, Vec U)
1093: {
1094: TS_RosW *ros = (TS_RosW *)ts->data;
1095: PetscInt s = ros->tableau->s, pinterp = ros->tableau->pinterp, i, j;
1096: PetscReal h;
1097: PetscReal tt, t;
1098: PetscScalar *bt;
1099: const PetscReal *Bt = ros->tableau->binterpt;
1100: const PetscReal *GammaInv = ros->tableau->GammaInv;
1101: PetscScalar *w = ros->work;
1102: Vec *Y = ros->Y;
1104: PetscFunctionBegin;
1105: PetscCheck(Bt, PetscObjectComm((PetscObject)ts), PETSC_ERR_SUP, "TSRosW %s does not have an interpolation formula", ros->tableau->name);
1107: switch (ros->status) {
1108: case TS_STEP_INCOMPLETE:
1109: case TS_STEP_PENDING:
1110: h = ts->time_step;
1111: t = (itime - ts->ptime) / h;
1112: break;
1113: case TS_STEP_COMPLETE:
1114: h = ts->ptime - ts->ptime_prev;
1115: t = (itime - ts->ptime) / h + 1; /* In the interval [0,1] */
1116: break;
1117: default:
1118: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_PLIB, "Invalid TSStepStatus");
1119: }
1120: PetscCall(PetscMalloc1(s, &bt));
1121: for (i = 0; i < s; i++) bt[i] = 0;
1122: for (j = 0, tt = t; j < pinterp; j++, tt *= t) {
1123: for (i = 0; i < s; i++) bt[i] += Bt[i * pinterp + j] * tt;
1124: }
1126: /* y(t+tt*h) = y(t) + Sum bt(tt) * GammaInv * Ydot */
1127: /* U <- 0*/
1128: PetscCall(VecZeroEntries(U));
1129: /* U <- Sum bt_i * GammaInv(i,1:i) * Y(1:i) */
1130: for (j = 0; j < s; j++) w[j] = 0;
1131: for (j = 0; j < s; j++) {
1132: for (i = j; i < s; i++) w[j] += bt[i] * GammaInv[i * s + j];
1133: }
1134: PetscCall(VecMAXPY(U, i, w, Y));
1135: /* U <- y(t) + U */
1136: PetscCall(VecAXPY(U, 1, ros->vec_sol_prev));
1138: PetscCall(PetscFree(bt));
1139: PetscFunctionReturn(PETSC_SUCCESS);
1140: }
1142: /*------------------------------------------------------------*/
1144: static PetscErrorCode TSRosWTableauReset(TS ts)
1145: {
1146: TS_RosW *ros = (TS_RosW *)ts->data;
1147: RosWTableau tab = ros->tableau;
1149: PetscFunctionBegin;
1150: if (!tab) PetscFunctionReturn(PETSC_SUCCESS);
1151: PetscCall(VecDestroyVecs(tab->s, &ros->Y));
1152: PetscCall(PetscFree(ros->work));
1153: PetscFunctionReturn(PETSC_SUCCESS);
1154: }
1156: static PetscErrorCode TSReset_RosW(TS ts)
1157: {
1158: TS_RosW *ros = (TS_RosW *)ts->data;
1160: PetscFunctionBegin;
1161: PetscCall(TSRosWTableauReset(ts));
1162: PetscCall(VecDestroy(&ros->Ydot));
1163: PetscCall(VecDestroy(&ros->Ystage));
1164: PetscCall(VecDestroy(&ros->Zdot));
1165: PetscCall(VecDestroy(&ros->Zstage));
1166: PetscCall(VecDestroy(&ros->vec_sol_prev));
1167: PetscFunctionReturn(PETSC_SUCCESS);
1168: }
1170: static PetscErrorCode TSRosWGetVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1171: {
1172: TS_RosW *rw = (TS_RosW *)ts->data;
1174: PetscFunctionBegin;
1175: if (Ydot) {
1176: if (dm && dm != ts->dm) {
1177: PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1178: } else *Ydot = rw->Ydot;
1179: }
1180: if (Zdot) {
1181: if (dm && dm != ts->dm) {
1182: PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1183: } else *Zdot = rw->Zdot;
1184: }
1185: if (Ystage) {
1186: if (dm && dm != ts->dm) {
1187: PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1188: } else *Ystage = rw->Ystage;
1189: }
1190: if (Zstage) {
1191: if (dm && dm != ts->dm) {
1192: PetscCall(DMGetNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1193: } else *Zstage = rw->Zstage;
1194: }
1195: PetscFunctionReturn(PETSC_SUCCESS);
1196: }
1198: static PetscErrorCode TSRosWRestoreVecs(TS ts, DM dm, Vec *Ydot, Vec *Zdot, Vec *Ystage, Vec *Zstage)
1199: {
1200: PetscFunctionBegin;
1201: if (Ydot) {
1202: if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ydot", Ydot));
1203: }
1204: if (Zdot) {
1205: if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zdot", Zdot));
1206: }
1207: if (Ystage) {
1208: if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Ystage", Ystage));
1209: }
1210: if (Zstage) {
1211: if (dm && dm != ts->dm) PetscCall(DMRestoreNamedGlobalVector(dm, "TSRosW_Zstage", Zstage));
1212: }
1213: PetscFunctionReturn(PETSC_SUCCESS);
1214: }
1216: static PetscErrorCode DMCoarsenHook_TSRosW(DM fine, DM coarse, void *ctx)
1217: {
1218: PetscFunctionBegin;
1219: PetscFunctionReturn(PETSC_SUCCESS);
1220: }
1222: static PetscErrorCode DMRestrictHook_TSRosW(DM fine, Mat restrct, Vec rscale, Mat inject, DM coarse, void *ctx)
1223: {
1224: TS ts = (TS)ctx;
1225: Vec Ydot, Zdot, Ystage, Zstage;
1226: Vec Ydotc, Zdotc, Ystagec, Zstagec;
1228: PetscFunctionBegin;
1229: PetscCall(TSRosWGetVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
1230: PetscCall(TSRosWGetVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
1231: PetscCall(MatRestrict(restrct, Ydot, Ydotc));
1232: PetscCall(VecPointwiseMult(Ydotc, rscale, Ydotc));
1233: PetscCall(MatRestrict(restrct, Ystage, Ystagec));
1234: PetscCall(VecPointwiseMult(Ystagec, rscale, Ystagec));
1235: PetscCall(MatRestrict(restrct, Zdot, Zdotc));
1236: PetscCall(VecPointwiseMult(Zdotc, rscale, Zdotc));
1237: PetscCall(MatRestrict(restrct, Zstage, Zstagec));
1238: PetscCall(VecPointwiseMult(Zstagec, rscale, Zstagec));
1239: PetscCall(TSRosWRestoreVecs(ts, fine, &Ydot, &Ystage, &Zdot, &Zstage));
1240: PetscCall(TSRosWRestoreVecs(ts, coarse, &Ydotc, &Ystagec, &Zdotc, &Zstagec));
1241: PetscFunctionReturn(PETSC_SUCCESS);
1242: }
1244: static PetscErrorCode DMSubDomainHook_TSRosW(DM fine, DM coarse, void *ctx)
1245: {
1246: PetscFunctionBegin;
1247: PetscFunctionReturn(PETSC_SUCCESS);
1248: }
1250: static PetscErrorCode DMSubDomainRestrictHook_TSRosW(DM dm, VecScatter gscat, VecScatter lscat, DM subdm, void *ctx)
1251: {
1252: TS ts = (TS)ctx;
1253: Vec Ydot, Zdot, Ystage, Zstage;
1254: Vec Ydots, Zdots, Ystages, Zstages;
1256: PetscFunctionBegin;
1257: PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
1258: PetscCall(TSRosWGetVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
1260: PetscCall(VecScatterBegin(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
1261: PetscCall(VecScatterEnd(gscat, Ydot, Ydots, INSERT_VALUES, SCATTER_FORWARD));
1263: PetscCall(VecScatterBegin(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
1264: PetscCall(VecScatterEnd(gscat, Ystage, Ystages, INSERT_VALUES, SCATTER_FORWARD));
1266: PetscCall(VecScatterBegin(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
1267: PetscCall(VecScatterEnd(gscat, Zdot, Zdots, INSERT_VALUES, SCATTER_FORWARD));
1269: PetscCall(VecScatterBegin(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
1270: PetscCall(VecScatterEnd(gscat, Zstage, Zstages, INSERT_VALUES, SCATTER_FORWARD));
1272: PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Ystage, &Zdot, &Zstage));
1273: PetscCall(TSRosWRestoreVecs(ts, subdm, &Ydots, &Ystages, &Zdots, &Zstages));
1274: PetscFunctionReturn(PETSC_SUCCESS);
1275: }
1277: /*
1278: This defines the nonlinear equation that is to be solved with SNES
1279: G(U) = F[t0+Theta*dt, U, (U-U0)*shift] = 0
1280: */
1281: static PetscErrorCode SNESTSFormFunction_RosW(SNES snes, Vec U, Vec F, TS ts)
1282: {
1283: TS_RosW *ros = (TS_RosW *)ts->data;
1284: Vec Ydot, Zdot, Ystage, Zstage;
1285: PetscReal shift = ros->scoeff / ts->time_step;
1286: DM dm, dmsave;
1288: PetscFunctionBegin;
1289: PetscCall(SNESGetDM(snes, &dm));
1290: PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1291: PetscCall(VecWAXPY(Ydot, shift, U, Zdot)); /* Ydot = shift*U + Zdot */
1292: PetscCall(VecWAXPY(Ystage, 1.0, U, Zstage)); /* Ystage = U + Zstage */
1293: dmsave = ts->dm;
1294: ts->dm = dm;
1295: PetscCall(TSComputeIFunction(ts, ros->stage_time, Ystage, Ydot, F, PETSC_FALSE));
1296: ts->dm = dmsave;
1297: PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1298: PetscFunctionReturn(PETSC_SUCCESS);
1299: }
1301: static PetscErrorCode SNESTSFormJacobian_RosW(SNES snes, Vec U, Mat A, Mat B, TS ts)
1302: {
1303: TS_RosW *ros = (TS_RosW *)ts->data;
1304: Vec Ydot, Zdot, Ystage, Zstage;
1305: PetscReal shift = ros->scoeff / ts->time_step;
1306: DM dm, dmsave;
1308: PetscFunctionBegin;
1309: /* ros->Ydot and ros->Ystage have already been computed in SNESTSFormFunction_RosW (SNES guarantees this) */
1310: PetscCall(SNESGetDM(snes, &dm));
1311: PetscCall(TSRosWGetVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1312: dmsave = ts->dm;
1313: ts->dm = dm;
1314: PetscCall(TSComputeIJacobian(ts, ros->stage_time, Ystage, Ydot, shift, A, B, PETSC_TRUE));
1315: ts->dm = dmsave;
1316: PetscCall(TSRosWRestoreVecs(ts, dm, &Ydot, &Zdot, &Ystage, &Zstage));
1317: PetscFunctionReturn(PETSC_SUCCESS);
1318: }
1320: static PetscErrorCode TSRosWTableauSetUp(TS ts)
1321: {
1322: TS_RosW *ros = (TS_RosW *)ts->data;
1323: RosWTableau tab = ros->tableau;
1325: PetscFunctionBegin;
1326: PetscCall(VecDuplicateVecs(ts->vec_sol, tab->s, &ros->Y));
1327: PetscCall(PetscMalloc1(tab->s, &ros->work));
1328: PetscFunctionReturn(PETSC_SUCCESS);
1329: }
1331: static PetscErrorCode TSSetUp_RosW(TS ts)
1332: {
1333: TS_RosW *ros = (TS_RosW *)ts->data;
1334: DM dm;
1335: SNES snes;
1336: TSRHSJacobian rhsjacobian;
1338: PetscFunctionBegin;
1339: PetscCall(TSRosWTableauSetUp(ts));
1340: PetscCall(VecDuplicate(ts->vec_sol, &ros->Ydot));
1341: PetscCall(VecDuplicate(ts->vec_sol, &ros->Ystage));
1342: PetscCall(VecDuplicate(ts->vec_sol, &ros->Zdot));
1343: PetscCall(VecDuplicate(ts->vec_sol, &ros->Zstage));
1344: PetscCall(VecDuplicate(ts->vec_sol, &ros->vec_sol_prev));
1345: PetscCall(TSGetDM(ts, &dm));
1346: PetscCall(DMCoarsenHookAdd(dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
1347: PetscCall(DMSubDomainHookAdd(dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1348: /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1349: PetscCall(TSGetSNES(ts, &snes));
1350: if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
1351: PetscCall(DMTSGetRHSJacobian(dm, &rhsjacobian, NULL));
1352: if (rhsjacobian == TSComputeRHSJacobianConstant) {
1353: Mat Amat, Pmat;
1355: /* Set the SNES matrix to be different from the RHS matrix because there is no way to reconstruct shift*M-J */
1356: PetscCall(SNESGetJacobian(snes, &Amat, &Pmat, NULL, NULL));
1357: if (Amat && Amat == ts->Arhs) {
1358: if (Amat == Pmat) {
1359: PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
1360: PetscCall(SNESSetJacobian(snes, Amat, Amat, NULL, NULL));
1361: } else {
1362: PetscCall(MatDuplicate(ts->Arhs, MAT_COPY_VALUES, &Amat));
1363: PetscCall(SNESSetJacobian(snes, Amat, NULL, NULL, NULL));
1364: if (Pmat && Pmat == ts->Brhs) {
1365: PetscCall(MatDuplicate(ts->Brhs, MAT_COPY_VALUES, &Pmat));
1366: PetscCall(SNESSetJacobian(snes, NULL, Pmat, NULL, NULL));
1367: PetscCall(MatDestroy(&Pmat));
1368: }
1369: }
1370: PetscCall(MatDestroy(&Amat));
1371: }
1372: }
1373: PetscFunctionReturn(PETSC_SUCCESS);
1374: }
1375: /*------------------------------------------------------------*/
1377: static PetscErrorCode TSSetFromOptions_RosW(TS ts, PetscOptionItems *PetscOptionsObject)
1378: {
1379: TS_RosW *ros = (TS_RosW *)ts->data;
1380: SNES snes;
1382: PetscFunctionBegin;
1383: PetscOptionsHeadBegin(PetscOptionsObject, "RosW ODE solver options");
1384: {
1385: RosWTableauLink link;
1386: PetscInt count, choice;
1387: PetscBool flg;
1388: const char **namelist;
1390: for (link = RosWTableauList, count = 0; link; link = link->next, count++)
1391: ;
1392: PetscCall(PetscMalloc1(count, (char ***)&namelist));
1393: for (link = RosWTableauList, count = 0; link; link = link->next, count++) namelist[count] = link->tab.name;
1394: PetscCall(PetscOptionsEList("-ts_rosw_type", "Family of Rosenbrock-W method", "TSRosWSetType", (const char *const *)namelist, count, ros->tableau->name, &choice, &flg));
1395: if (flg) PetscCall(TSRosWSetType(ts, namelist[choice]));
1396: PetscCall(PetscFree(namelist));
1398: PetscCall(PetscOptionsBool("-ts_rosw_recompute_jacobian", "Recompute the Jacobian at each stage", "TSRosWSetRecomputeJacobian", ros->recompute_jacobian, &ros->recompute_jacobian, NULL));
1399: }
1400: PetscOptionsHeadEnd();
1401: /* Rosenbrock methods are linearly implicit, so set that unless the user has specifically asked for something else */
1402: PetscCall(TSGetSNES(ts, &snes));
1403: if (!((PetscObject)snes)->type_name) PetscCall(SNESSetType(snes, SNESKSPONLY));
1404: PetscFunctionReturn(PETSC_SUCCESS);
1405: }
1407: static PetscErrorCode TSView_RosW(TS ts, PetscViewer viewer)
1408: {
1409: TS_RosW *ros = (TS_RosW *)ts->data;
1410: PetscBool iascii;
1412: PetscFunctionBegin;
1413: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
1414: if (iascii) {
1415: RosWTableau tab = ros->tableau;
1416: TSRosWType rostype;
1417: char buf[512];
1418: PetscInt i;
1419: PetscReal abscissa[512];
1420: PetscCall(TSRosWGetType(ts, &rostype));
1421: PetscCall(PetscViewerASCIIPrintf(viewer, " Rosenbrock-W %s\n", rostype));
1422: PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, tab->ASum));
1423: PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa of A = %s\n", buf));
1424: for (i = 0; i < tab->s; i++) abscissa[i] = tab->ASum[i] + tab->Gamma[i];
1425: PetscCall(PetscFormatRealArray(buf, sizeof(buf), "% 8.6f", tab->s, abscissa));
1426: PetscCall(PetscViewerASCIIPrintf(viewer, " Abscissa of A+Gamma = %s\n", buf));
1427: }
1428: PetscFunctionReturn(PETSC_SUCCESS);
1429: }
1431: static PetscErrorCode TSLoad_RosW(TS ts, PetscViewer viewer)
1432: {
1433: SNES snes;
1434: TSAdapt adapt;
1436: PetscFunctionBegin;
1437: PetscCall(TSGetAdapt(ts, &adapt));
1438: PetscCall(TSAdaptLoad(adapt, viewer));
1439: PetscCall(TSGetSNES(ts, &snes));
1440: PetscCall(SNESLoad(snes, viewer));
1441: /* function and Jacobian context for SNES when used with TS is always ts object */
1442: PetscCall(SNESSetFunction(snes, NULL, NULL, ts));
1443: PetscCall(SNESSetJacobian(snes, NULL, NULL, NULL, ts));
1444: PetscFunctionReturn(PETSC_SUCCESS);
1445: }
1447: /*@C
1448: TSRosWSetType - Set the type of Rosenbrock-W, `TSROSW`, scheme
1450: Logically collective
1452: Input Parameters:
1453: + ts - timestepping context
1454: - roswtype - type of Rosenbrock-W scheme
1456: Level: beginner
1458: .seealso: [](chapter_ts), `TSRosWGetType()`, `TSROSW`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`, `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWARK3`
1459: @*/
1460: PetscErrorCode TSRosWSetType(TS ts, TSRosWType roswtype)
1461: {
1462: PetscFunctionBegin;
1465: PetscTryMethod(ts, "TSRosWSetType_C", (TS, TSRosWType), (ts, roswtype));
1466: PetscFunctionReturn(PETSC_SUCCESS);
1467: }
1469: /*@C
1470: TSRosWGetType - Get the type of Rosenbrock-W scheme
1472: Logically collective
1474: Input Parameter:
1475: . ts - timestepping context
1477: Output Parameter:
1478: . rostype - type of Rosenbrock-W scheme
1480: Level: intermediate
1482: .seealso: [](chapter_ts), `TSRosWType`, `TSRosWSetType()`
1483: @*/
1484: PetscErrorCode TSRosWGetType(TS ts, TSRosWType *rostype)
1485: {
1486: PetscFunctionBegin;
1488: PetscUseMethod(ts, "TSRosWGetType_C", (TS, TSRosWType *), (ts, rostype));
1489: PetscFunctionReturn(PETSC_SUCCESS);
1490: }
1492: /*@C
1493: TSRosWSetRecomputeJacobian - Set whether to recompute the Jacobian at each stage. The default is to update the Jacobian once per step.
1495: Logically collective
1497: Input Parameters:
1498: + ts - timestepping context
1499: - flg - `PETSC_TRUE` to recompute the Jacobian at each stage
1501: Level: intermediate
1503: .seealso: [](chapter_ts), `TSRosWType`, `TSRosWGetType()`
1504: @*/
1505: PetscErrorCode TSRosWSetRecomputeJacobian(TS ts, PetscBool flg)
1506: {
1507: PetscFunctionBegin;
1509: PetscTryMethod(ts, "TSRosWSetRecomputeJacobian_C", (TS, PetscBool), (ts, flg));
1510: PetscFunctionReturn(PETSC_SUCCESS);
1511: }
1513: static PetscErrorCode TSRosWGetType_RosW(TS ts, TSRosWType *rostype)
1514: {
1515: TS_RosW *ros = (TS_RosW *)ts->data;
1517: PetscFunctionBegin;
1518: *rostype = ros->tableau->name;
1519: PetscFunctionReturn(PETSC_SUCCESS);
1520: }
1522: static PetscErrorCode TSRosWSetType_RosW(TS ts, TSRosWType rostype)
1523: {
1524: TS_RosW *ros = (TS_RosW *)ts->data;
1525: PetscBool match;
1526: RosWTableauLink link;
1528: PetscFunctionBegin;
1529: if (ros->tableau) {
1530: PetscCall(PetscStrcmp(ros->tableau->name, rostype, &match));
1531: if (match) PetscFunctionReturn(PETSC_SUCCESS);
1532: }
1533: for (link = RosWTableauList; link; link = link->next) {
1534: PetscCall(PetscStrcmp(link->tab.name, rostype, &match));
1535: if (match) {
1536: if (ts->setupcalled) PetscCall(TSRosWTableauReset(ts));
1537: ros->tableau = &link->tab;
1538: if (ts->setupcalled) PetscCall(TSRosWTableauSetUp(ts));
1539: ts->default_adapt_type = ros->tableau->bembed ? TSADAPTBASIC : TSADAPTNONE;
1540: PetscFunctionReturn(PETSC_SUCCESS);
1541: }
1542: }
1543: SETERRQ(PetscObjectComm((PetscObject)ts), PETSC_ERR_ARG_UNKNOWN_TYPE, "Could not find '%s'", rostype);
1544: }
1546: static PetscErrorCode TSRosWSetRecomputeJacobian_RosW(TS ts, PetscBool flg)
1547: {
1548: TS_RosW *ros = (TS_RosW *)ts->data;
1550: PetscFunctionBegin;
1551: ros->recompute_jacobian = flg;
1552: PetscFunctionReturn(PETSC_SUCCESS);
1553: }
1555: static PetscErrorCode TSDestroy_RosW(TS ts)
1556: {
1557: PetscFunctionBegin;
1558: PetscCall(TSReset_RosW(ts));
1559: if (ts->dm) {
1560: PetscCall(DMCoarsenHookRemove(ts->dm, DMCoarsenHook_TSRosW, DMRestrictHook_TSRosW, ts));
1561: PetscCall(DMSubDomainHookRemove(ts->dm, DMSubDomainHook_TSRosW, DMSubDomainRestrictHook_TSRosW, ts));
1562: }
1563: PetscCall(PetscFree(ts->data));
1564: PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", NULL));
1565: PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", NULL));
1566: PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", NULL));
1567: PetscFunctionReturn(PETSC_SUCCESS);
1568: }
1570: /* ------------------------------------------------------------ */
1571: /*MC
1572: TSROSW - ODE solver using Rosenbrock-W schemes
1574: These methods are intended for problems with well-separated time scales, especially when a slow scale is strongly
1575: nonlinear such that it is expensive to solve with a fully implicit method. The user should provide the stiff part
1576: of the equation using `TSSetIFunction()` and the non-stiff part with `TSSetRHSFunction()`.
1578: Level: beginner
1580: Notes:
1581: This method currently only works with autonomous ODE and DAE.
1583: Consider trying `TSARKIMEX` if the stiff part is strongly nonlinear.
1585: 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
1587: Developer Notes:
1588: Rosenbrock-W methods are typically specified for autonomous ODE
1590: $ udot = f(u)
1592: by the stage equations
1594: $ k_i = h f(u_0 + sum_j alpha_ij k_j) + h J sum_j gamma_ij k_j
1596: and step completion formula
1598: $ u_1 = u_0 + sum_j b_j k_j
1600: with step size h and coefficients alpha_ij, gamma_ij, and b_i. Implementing the method in this form would require f(u)
1601: and the Jacobian J to be available, in addition to the shifted matrix I - h gamma_ii J. Following Hairer and Wanner,
1602: we define new variables for the stage equations
1604: $ y_i = gamma_ij k_j
1606: The k_j can be recovered because Gamma is invertible. Let C be the lower triangular part of Gamma^{-1} and define
1608: $ A = Alpha Gamma^{-1}, bt^T = b^T Gamma^{-1}
1610: to rewrite the method as
1612: $ [M/(h gamma_ii) - J] y_i = f(u_0 + sum_j a_ij y_j) + M sum_j (c_ij/h) y_j
1613: $ u_1 = u_0 + sum_j bt_j y_j
1615: where we have introduced the mass matrix M. Continue by defining
1617: $ ydot_i = 1/(h gamma_ii) y_i - sum_j (c_ij/h) y_j
1619: or, more compactly in tensor notation
1621: $ Ydot = 1/h (Gamma^{-1} \otimes I) Y .
1623: Note that Gamma^{-1} is lower triangular. With this definition of Ydot in terms of known quantities and the current
1624: stage y_i, the stage equations reduce to performing one Newton step (typically with a lagged Jacobian) on the
1625: equation
1627: $ g(u_0 + sum_j a_ij y_j + y_i, ydot_i) = 0
1629: with initial guess y_i = 0.
1631: .seealso: [](chapter_ts), `TSCreate()`, `TS`, `TSSetType()`, `TSRosWSetType()`, `TSRosWRegister()`, `TSROSWTHETA1`, `TSROSWTHETA2`, `TSROSW2M`, `TSROSW2P`, `TSROSWRA3PW`, `TSROSWRA34PW2`, `TSROSWRODAS3`,
1632: `TSROSWSANDU3`, `TSROSWASSP3P3S1C`, `TSROSWLASSP3P4S2C`, `TSROSWLLSSP3P4S2C`, `TSROSWGRK4T`, `TSROSWSHAMP4`, `TSROSWVELDD4`, `TSROSW4L`, `TSType`
1633: M*/
1634: PETSC_EXTERN PetscErrorCode TSCreate_RosW(TS ts)
1635: {
1636: TS_RosW *ros;
1638: PetscFunctionBegin;
1639: PetscCall(TSRosWInitializePackage());
1641: ts->ops->reset = TSReset_RosW;
1642: ts->ops->destroy = TSDestroy_RosW;
1643: ts->ops->view = TSView_RosW;
1644: ts->ops->load = TSLoad_RosW;
1645: ts->ops->setup = TSSetUp_RosW;
1646: ts->ops->step = TSStep_RosW;
1647: ts->ops->interpolate = TSInterpolate_RosW;
1648: ts->ops->evaluatestep = TSEvaluateStep_RosW;
1649: ts->ops->rollback = TSRollBack_RosW;
1650: ts->ops->setfromoptions = TSSetFromOptions_RosW;
1651: ts->ops->snesfunction = SNESTSFormFunction_RosW;
1652: ts->ops->snesjacobian = SNESTSFormJacobian_RosW;
1654: ts->usessnes = PETSC_TRUE;
1656: PetscCall(PetscNew(&ros));
1657: ts->data = (void *)ros;
1659: PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWGetType_C", TSRosWGetType_RosW));
1660: PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetType_C", TSRosWSetType_RosW));
1661: PetscCall(PetscObjectComposeFunction((PetscObject)ts, "TSRosWSetRecomputeJacobian_C", TSRosWSetRecomputeJacobian_RosW));
1663: PetscCall(TSRosWSetType(ts, TSRosWDefault));
1664: PetscFunctionReturn(PETSC_SUCCESS);
1665: }