xref: /petsc/src/ts/tutorials/ex31.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1c4762a1bSJed Brown static char help[] = "Solves the ordinary differential equations (IVPs) using explicit and implicit time-integration methods.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown 
5c4762a1bSJed Brown   Useful command line parameters:
6c4762a1bSJed Brown   -problem <hull1972a1>: choose which problem to solve (see references
7c4762a1bSJed Brown                       for complete listing of problems).
8c4762a1bSJed Brown   -ts_type <euler>: specify time-integrator
9c4762a1bSJed Brown   -ts_adapt_type <basic>: specify time-step adapting (none,basic,advanced)
10c4762a1bSJed Brown   -refinement_levels <1>: number of refinement levels for convergence analysis
11c4762a1bSJed Brown   -refinement_factor <2.0>: factor to refine time step size by for convergence analysis
12c4762a1bSJed Brown   -dt <0.01>: specify time step (initial time step for convergence analysis)
13c4762a1bSJed Brown 
14c4762a1bSJed Brown */
15c4762a1bSJed Brown 
16c4762a1bSJed Brown /*
17c4762a1bSJed Brown List of cases and their names in the code:-
18c4762a1bSJed Brown   From Hull, T.E., Enright, W.H., Fellen, B.M., and Sedgwick, A.E.,
19c4762a1bSJed Brown       "Comparing Numerical Methods for Ordinary Differential
20c4762a1bSJed Brown        Equations", SIAM J. Numer. Anal., 9(4), 1972, pp. 603 - 635
21c4762a1bSJed Brown     A1 -> "hull1972a1" (exact solution available)
22c4762a1bSJed Brown     A2 -> "hull1972a2" (exact solution available)
23c4762a1bSJed Brown     A3 -> "hull1972a3" (exact solution available)
24c4762a1bSJed Brown     A4 -> "hull1972a4" (exact solution available)
25c4762a1bSJed Brown     A5 -> "hull1972a5"
26c4762a1bSJed Brown     B1 -> "hull1972b1"
27c4762a1bSJed Brown     B2 -> "hull1972b2"
28c4762a1bSJed Brown     B3 -> "hull1972b3"
29c4762a1bSJed Brown     B4 -> "hull1972b4"
30c4762a1bSJed Brown     B5 -> "hull1972b5"
31c4762a1bSJed Brown     C1 -> "hull1972c1"
32c4762a1bSJed Brown     C2 -> "hull1972c2"
33c4762a1bSJed Brown     C3 -> "hull1972c3"
34c4762a1bSJed Brown     C4 -> "hull1972c4"
35c4762a1bSJed Brown 
36c4762a1bSJed Brown  From Constantinescu, E. "Estimating Global Errors in Time Stepping" ArXiv e-prints,
37c4762a1bSJed Brown        https://arxiv.org/abs/1503.05166, 2016
38c4762a1bSJed Brown 
39c4762a1bSJed Brown     Kulikov2013I -> "kulik2013i"
40c4762a1bSJed Brown 
41c4762a1bSJed Brown */
42c4762a1bSJed Brown 
43c4762a1bSJed Brown #include <petscts.h>
44c4762a1bSJed Brown 
45c4762a1bSJed Brown /* Function declarations */
46c4762a1bSJed Brown PetscErrorCode (*RHSFunction)(TS, PetscReal, Vec, Vec, void *);
47c4762a1bSJed Brown PetscErrorCode (*RHSJacobian)(TS, PetscReal, Vec, Mat, Mat, void *);
48c4762a1bSJed Brown PetscErrorCode (*IFunction)(TS, PetscReal, Vec, Vec, Vec, void *);
49c4762a1bSJed Brown PetscErrorCode (*IJacobian)(TS, PetscReal, Vec, Vec, PetscReal, Mat, Mat, void *);
50c4762a1bSJed Brown 
51c4762a1bSJed Brown /* Returns the size of the system of equations depending on problem specification */
52*9371c9d4SSatish Balay PetscInt GetSize(const char *p) {
537510d9b0SBarry Smith   PetscFunctionBeginUser;
54*9371c9d4SSatish Balay   if ((!strcmp(p, "hull1972a1")) || (!strcmp(p, "hull1972a2")) || (!strcmp(p, "hull1972a3")) || (!strcmp(p, "hull1972a4")) || (!strcmp(p, "hull1972a5"))) PetscFunctionReturn(1);
55c4762a1bSJed Brown   else if (!strcmp(p, "hull1972b1")) PetscFunctionReturn(2);
56*9371c9d4SSatish Balay   else if ((!strcmp(p, "hull1972b2")) || (!strcmp(p, "hull1972b3")) || (!strcmp(p, "hull1972b4")) || (!strcmp(p, "hull1972b5"))) PetscFunctionReturn(3);
57c4762a1bSJed Brown   else if ((!strcmp(p, "kulik2013i"))) PetscFunctionReturn(4);
58*9371c9d4SSatish Balay   else if ((!strcmp(p, "hull1972c1")) || (!strcmp(p, "hull1972c2")) || (!strcmp(p, "hull1972c3"))) PetscFunctionReturn(10);
59c4762a1bSJed Brown   else if (!strcmp(p, "hull1972c4")) PetscFunctionReturn(51);
60c4762a1bSJed Brown   else PetscFunctionReturn(-1);
61c4762a1bSJed Brown }
62c4762a1bSJed Brown 
63c4762a1bSJed Brown /****************************************************************/
64c4762a1bSJed Brown 
65c4762a1bSJed Brown /* Problem specific functions */
66c4762a1bSJed Brown 
67c4762a1bSJed Brown /* Hull, 1972, Problem A1 */
68c4762a1bSJed Brown 
69*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
70c4762a1bSJed Brown   PetscScalar       *f;
71c4762a1bSJed Brown   const PetscScalar *y;
72c4762a1bSJed Brown 
737510d9b0SBarry Smith   PetscFunctionBeginUser;
749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
76c4762a1bSJed Brown   f[0] = -y[0];
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
789566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
79c4762a1bSJed Brown   PetscFunctionReturn(0);
80c4762a1bSJed Brown }
81c4762a1bSJed Brown 
82*9371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) {
83c4762a1bSJed Brown   const PetscScalar *y;
84c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
85c4762a1bSJed Brown   PetscScalar        value = -1.0;
86c4762a1bSJed Brown 
877510d9b0SBarry Smith   PetscFunctionBeginUser;
889566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
899566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
97c4762a1bSJed Brown   const PetscScalar *y;
98c4762a1bSJed Brown   PetscScalar       *f;
99c4762a1bSJed Brown 
1007510d9b0SBarry Smith   PetscFunctionBeginUser;
1019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1029566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
103c4762a1bSJed Brown   f[0] = -y[0];
1049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
1059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
106c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
1079566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
108c4762a1bSJed Brown   PetscFunctionReturn(0);
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
111*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
112c4762a1bSJed Brown   const PetscScalar *y;
113c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
114c4762a1bSJed Brown   PetscScalar        value = a - 1.0;
115c4762a1bSJed Brown 
1167510d9b0SBarry Smith   PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1189566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
1199566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1209566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
122c4762a1bSJed Brown   PetscFunctionReturn(0);
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown /* Hull, 1972, Problem A2 */
126c4762a1bSJed Brown 
127*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
128c4762a1bSJed Brown   const PetscScalar *y;
129c4762a1bSJed Brown   PetscScalar       *f;
130c4762a1bSJed Brown 
1317510d9b0SBarry Smith   PetscFunctionBeginUser;
1329566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1339566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
134c4762a1bSJed Brown   f[0] = -0.5 * y[0] * y[0] * y[0];
1359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
1369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
137c4762a1bSJed Brown   PetscFunctionReturn(0);
138c4762a1bSJed Brown }
139c4762a1bSJed Brown 
140*9371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) {
141c4762a1bSJed Brown   const PetscScalar *y;
142c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
143c4762a1bSJed Brown   PetscScalar        value;
144c4762a1bSJed Brown 
1457510d9b0SBarry Smith   PetscFunctionBeginUser;
1469566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
147c4762a1bSJed Brown   value = -0.5 * 3.0 * y[0] * y[0];
1489566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
1499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
152c4762a1bSJed Brown   PetscFunctionReturn(0);
153c4762a1bSJed Brown }
154c4762a1bSJed Brown 
155*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
156c4762a1bSJed Brown   PetscScalar       *f;
157c4762a1bSJed Brown   const PetscScalar *y;
158c4762a1bSJed Brown 
1597510d9b0SBarry Smith   PetscFunctionBeginUser;
1609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
162c4762a1bSJed Brown   f[0] = -0.5 * y[0] * y[0] * y[0];
1639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
1649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
165c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
1669566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
167c4762a1bSJed Brown   PetscFunctionReturn(0);
168c4762a1bSJed Brown }
169c4762a1bSJed Brown 
170*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
171c4762a1bSJed Brown   const PetscScalar *y;
172c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
173c4762a1bSJed Brown   PetscScalar        value;
174c4762a1bSJed Brown 
1757510d9b0SBarry Smith   PetscFunctionBeginUser;
1769566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
177c4762a1bSJed Brown   value = a + 0.5 * 3.0 * y[0] * y[0];
1789566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
1799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1809566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
182c4762a1bSJed Brown   PetscFunctionReturn(0);
183c4762a1bSJed Brown }
184c4762a1bSJed Brown 
185c4762a1bSJed Brown /* Hull, 1972, Problem A3 */
186c4762a1bSJed Brown 
187*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
188c4762a1bSJed Brown   const PetscScalar *y;
189c4762a1bSJed Brown   PetscScalar       *f;
190c4762a1bSJed Brown 
1917510d9b0SBarry Smith   PetscFunctionBeginUser;
1929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
194c4762a1bSJed Brown   f[0] = y[0] * PetscCosReal(t);
1959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
1969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
197c4762a1bSJed Brown   PetscFunctionReturn(0);
198c4762a1bSJed Brown }
199c4762a1bSJed Brown 
200*9371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) {
201c4762a1bSJed Brown   const PetscScalar *y;
202c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
203c4762a1bSJed Brown   PetscScalar        value = PetscCosReal(t);
204c4762a1bSJed Brown 
2057510d9b0SBarry Smith   PetscFunctionBeginUser;
2069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2079566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
2089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2099566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
211c4762a1bSJed Brown   PetscFunctionReturn(0);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown 
214*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
215c4762a1bSJed Brown   PetscScalar       *f;
216c4762a1bSJed Brown   const PetscScalar *y;
217c4762a1bSJed Brown 
2187510d9b0SBarry Smith   PetscFunctionBeginUser;
2199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2209566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
221c4762a1bSJed Brown   f[0] = y[0] * PetscCosReal(t);
2229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
2239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
224c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
2259566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
226c4762a1bSJed Brown   PetscFunctionReturn(0);
227c4762a1bSJed Brown }
228c4762a1bSJed Brown 
229*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
230c4762a1bSJed Brown   const PetscScalar *y;
231c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
232c4762a1bSJed Brown   PetscScalar        value = a - PetscCosReal(t);
233c4762a1bSJed Brown 
2347510d9b0SBarry Smith   PetscFunctionBeginUser;
2359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2369566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
2379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
240c4762a1bSJed Brown   PetscFunctionReturn(0);
241c4762a1bSJed Brown }
242c4762a1bSJed Brown 
243c4762a1bSJed Brown /* Hull, 1972, Problem A4 */
244c4762a1bSJed Brown 
245*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
246c4762a1bSJed Brown   PetscScalar       *f;
247c4762a1bSJed Brown   const PetscScalar *y;
248c4762a1bSJed Brown 
2497510d9b0SBarry Smith   PetscFunctionBeginUser;
2509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
252c4762a1bSJed Brown   f[0] = (0.25 * y[0]) * (1.0 - 0.05 * y[0]);
2539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
2549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
255c4762a1bSJed Brown   PetscFunctionReturn(0);
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258*9371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) {
259c4762a1bSJed Brown   const PetscScalar *y;
260c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
261c4762a1bSJed Brown   PetscScalar        value;
262c4762a1bSJed Brown 
2637510d9b0SBarry Smith   PetscFunctionBeginUser;
2649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
265c4762a1bSJed Brown   value = 0.25 * (1.0 - 0.05 * y[0]) - (0.25 * y[0]) * 0.05;
2669566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
2679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
270c4762a1bSJed Brown   PetscFunctionReturn(0);
271c4762a1bSJed Brown }
272c4762a1bSJed Brown 
273*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
274c4762a1bSJed Brown   PetscScalar       *f;
275c4762a1bSJed Brown   const PetscScalar *y;
276c4762a1bSJed Brown 
2777510d9b0SBarry Smith   PetscFunctionBeginUser;
2789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
280c4762a1bSJed Brown   f[0] = (0.25 * y[0]) * (1.0 - 0.05 * y[0]);
2819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
2829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
283c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
2849566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
285c4762a1bSJed Brown   PetscFunctionReturn(0);
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
288*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
289c4762a1bSJed Brown   const PetscScalar *y;
290c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
291c4762a1bSJed Brown   PetscScalar        value;
292c4762a1bSJed Brown 
2937510d9b0SBarry Smith   PetscFunctionBeginUser;
2949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
295c4762a1bSJed Brown   value = a - 0.25 * (1.0 - 0.05 * y[0]) + (0.25 * y[0]) * 0.05;
2969566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
2979566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2989566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
300c4762a1bSJed Brown   PetscFunctionReturn(0);
301c4762a1bSJed Brown }
302c4762a1bSJed Brown 
303c4762a1bSJed Brown /* Hull, 1972, Problem A5 */
304c4762a1bSJed Brown 
305*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
306c4762a1bSJed Brown   PetscScalar       *f;
307c4762a1bSJed Brown   const PetscScalar *y;
308c4762a1bSJed Brown 
3097510d9b0SBarry Smith   PetscFunctionBeginUser;
3109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
3119566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
312c4762a1bSJed Brown   f[0] = (y[0] - t) / (y[0] + t);
3139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
3149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
315c4762a1bSJed Brown   PetscFunctionReturn(0);
316c4762a1bSJed Brown }
317c4762a1bSJed Brown 
318*9371c9d4SSatish Balay PetscErrorCode RHSJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) {
319c4762a1bSJed Brown   const PetscScalar *y;
320c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
321c4762a1bSJed Brown   PetscScalar        value;
322c4762a1bSJed Brown 
3237510d9b0SBarry Smith   PetscFunctionBeginUser;
3249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
325c4762a1bSJed Brown   value = 2 * t / ((t + y[0]) * (t + y[0]));
3269566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
3279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
330c4762a1bSJed Brown   PetscFunctionReturn(0);
331c4762a1bSJed Brown }
332c4762a1bSJed Brown 
333*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
334c4762a1bSJed Brown   PetscScalar       *f;
335c4762a1bSJed Brown   const PetscScalar *y;
336c4762a1bSJed Brown 
3377510d9b0SBarry Smith   PetscFunctionBeginUser;
3389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
3399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
340c4762a1bSJed Brown   f[0] = (y[0] - t) / (y[0] + t);
3419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
3429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
343c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
3449566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
345c4762a1bSJed Brown   PetscFunctionReturn(0);
346c4762a1bSJed Brown }
347c4762a1bSJed Brown 
348*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
349c4762a1bSJed Brown   const PetscScalar *y;
350c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
351c4762a1bSJed Brown   PetscScalar        value;
352c4762a1bSJed Brown 
3537510d9b0SBarry Smith   PetscFunctionBeginUser;
3549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
355c4762a1bSJed Brown   value = a - 2 * t / ((t + y[0]) * (t + y[0]));
3569566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
3579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
360c4762a1bSJed Brown   PetscFunctionReturn(0);
361c4762a1bSJed Brown }
362c4762a1bSJed Brown 
363c4762a1bSJed Brown /* Hull, 1972, Problem B1 */
364c4762a1bSJed Brown 
365*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
366c4762a1bSJed Brown   PetscScalar       *f;
367c4762a1bSJed Brown   const PetscScalar *y;
368c4762a1bSJed Brown 
3697510d9b0SBarry Smith   PetscFunctionBeginUser;
3709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
3719566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
372c4762a1bSJed Brown   f[0] = 2.0 * (y[0] - y[0] * y[1]);
373c4762a1bSJed Brown   f[1] = -(y[1] - y[0] * y[1]);
3749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
3759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
376c4762a1bSJed Brown   PetscFunctionReturn(0);
377c4762a1bSJed Brown }
378c4762a1bSJed Brown 
379*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
380c4762a1bSJed Brown   PetscScalar       *f;
381c4762a1bSJed Brown   const PetscScalar *y;
382c4762a1bSJed Brown 
3837510d9b0SBarry Smith   PetscFunctionBeginUser;
3849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
3859566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
386c4762a1bSJed Brown   f[0] = 2.0 * (y[0] - y[0] * y[1]);
387c4762a1bSJed Brown   f[1] = -(y[1] - y[0] * y[1]);
3889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
3899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
390c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
3919566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
392c4762a1bSJed Brown   PetscFunctionReturn(0);
393c4762a1bSJed Brown }
394c4762a1bSJed Brown 
395*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
396c4762a1bSJed Brown   const PetscScalar *y;
397c4762a1bSJed Brown   PetscInt           row[2] = {0, 1};
398c4762a1bSJed Brown   PetscScalar        value[2][2];
399c4762a1bSJed Brown 
4007510d9b0SBarry Smith   PetscFunctionBeginUser;
4019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
402*9371c9d4SSatish Balay   value[0][0] = a - 2.0 * (1.0 - y[1]);
403*9371c9d4SSatish Balay   value[0][1] = 2.0 * y[0];
404*9371c9d4SSatish Balay   value[1][0] = -y[1];
405*9371c9d4SSatish Balay   value[1][1] = a + 1.0 - y[0];
4069566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, &row[0], 2, &row[0], &value[0][0], INSERT_VALUES));
4079566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4089566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
410c4762a1bSJed Brown   PetscFunctionReturn(0);
411c4762a1bSJed Brown }
412c4762a1bSJed Brown 
413c4762a1bSJed Brown /* Hull, 1972, Problem B2 */
414c4762a1bSJed Brown 
415*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
416c4762a1bSJed Brown   PetscScalar       *f;
417c4762a1bSJed Brown   const PetscScalar *y;
418c4762a1bSJed Brown 
4197510d9b0SBarry Smith   PetscFunctionBeginUser;
4209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
422c4762a1bSJed Brown   f[0] = -y[0] + y[1];
423c4762a1bSJed Brown   f[1] = y[0] - 2.0 * y[1] + y[2];
424c4762a1bSJed Brown   f[2] = y[1] - y[2];
4259566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
4269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
427c4762a1bSJed Brown   PetscFunctionReturn(0);
428c4762a1bSJed Brown }
429c4762a1bSJed Brown 
430*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
431c4762a1bSJed Brown   PetscScalar       *f;
432c4762a1bSJed Brown   const PetscScalar *y;
433c4762a1bSJed Brown 
4347510d9b0SBarry Smith   PetscFunctionBeginUser;
4359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4369566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
437c4762a1bSJed Brown   f[0] = -y[0] + y[1];
438c4762a1bSJed Brown   f[1] = y[0] - 2.0 * y[1] + y[2];
439c4762a1bSJed Brown   f[2] = y[1] - y[2];
4409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
4419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
442c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
4439566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
444c4762a1bSJed Brown   PetscFunctionReturn(0);
445c4762a1bSJed Brown }
446c4762a1bSJed Brown 
447*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
448c4762a1bSJed Brown   const PetscScalar *y;
449c4762a1bSJed Brown   PetscInt           row[3] = {0, 1, 2};
450c4762a1bSJed Brown   PetscScalar        value[3][3];
451c4762a1bSJed Brown 
4527510d9b0SBarry Smith   PetscFunctionBeginUser;
4539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
454*9371c9d4SSatish Balay   value[0][0] = a + 1.0;
455*9371c9d4SSatish Balay   value[0][1] = -1.0;
456*9371c9d4SSatish Balay   value[0][2] = 0;
457*9371c9d4SSatish Balay   value[1][0] = -1.0;
458*9371c9d4SSatish Balay   value[1][1] = a + 2.0;
459*9371c9d4SSatish Balay   value[1][2] = -1.0;
460*9371c9d4SSatish Balay   value[2][0] = 0;
461*9371c9d4SSatish Balay   value[2][1] = -1.0;
462*9371c9d4SSatish Balay   value[2][2] = a + 1.0;
4639566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES));
4649566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
467c4762a1bSJed Brown   PetscFunctionReturn(0);
468c4762a1bSJed Brown }
469c4762a1bSJed Brown 
470c4762a1bSJed Brown /* Hull, 1972, Problem B3 */
471c4762a1bSJed Brown 
472*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
473c4762a1bSJed Brown   PetscScalar       *f;
474c4762a1bSJed Brown   const PetscScalar *y;
475c4762a1bSJed Brown 
4767510d9b0SBarry Smith   PetscFunctionBeginUser;
4779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4789566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
479c4762a1bSJed Brown   f[0] = -y[0];
480c4762a1bSJed Brown   f[1] = y[0] - y[1] * y[1];
481c4762a1bSJed Brown   f[2] = y[1] * y[1];
4829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
4839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
484c4762a1bSJed Brown   PetscFunctionReturn(0);
485c4762a1bSJed Brown }
486c4762a1bSJed Brown 
487*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
488c4762a1bSJed Brown   PetscScalar       *f;
489c4762a1bSJed Brown   const PetscScalar *y;
490c4762a1bSJed Brown 
4917510d9b0SBarry Smith   PetscFunctionBeginUser;
4929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
494c4762a1bSJed Brown   f[0] = -y[0];
495c4762a1bSJed Brown   f[1] = y[0] - y[1] * y[1];
496c4762a1bSJed Brown   f[2] = y[1] * y[1];
4979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
4989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
499c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
5009566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
501c4762a1bSJed Brown   PetscFunctionReturn(0);
502c4762a1bSJed Brown }
503c4762a1bSJed Brown 
504*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
505c4762a1bSJed Brown   const PetscScalar *y;
506c4762a1bSJed Brown   PetscInt           row[3] = {0, 1, 2};
507c4762a1bSJed Brown   PetscScalar        value[3][3];
508c4762a1bSJed Brown 
5097510d9b0SBarry Smith   PetscFunctionBeginUser;
5109566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
511*9371c9d4SSatish Balay   value[0][0] = a + 1.0;
512*9371c9d4SSatish Balay   value[0][1] = 0;
513*9371c9d4SSatish Balay   value[0][2] = 0;
514*9371c9d4SSatish Balay   value[1][0] = -1.0;
515*9371c9d4SSatish Balay   value[1][1] = a + 2.0 * y[1];
516*9371c9d4SSatish Balay   value[1][2] = 0;
517*9371c9d4SSatish Balay   value[2][0] = 0;
518*9371c9d4SSatish Balay   value[2][1] = -2.0 * y[1];
519*9371c9d4SSatish Balay   value[2][2] = a;
5209566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES));
5219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
524c4762a1bSJed Brown   PetscFunctionReturn(0);
525c4762a1bSJed Brown }
526c4762a1bSJed Brown 
527c4762a1bSJed Brown /* Hull, 1972, Problem B4 */
528c4762a1bSJed Brown 
529*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
530c4762a1bSJed Brown   PetscScalar       *f;
531c4762a1bSJed Brown   const PetscScalar *y;
532c4762a1bSJed Brown 
5337510d9b0SBarry Smith   PetscFunctionBeginUser;
5349566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
5359566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
536c4762a1bSJed Brown   f[0] = -y[1] - y[0] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
537c4762a1bSJed Brown   f[1] = y[0] - y[1] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
538c4762a1bSJed Brown   f[2] = y[0] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
5399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
5409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
541c4762a1bSJed Brown   PetscFunctionReturn(0);
542c4762a1bSJed Brown }
543c4762a1bSJed Brown 
544*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
545c4762a1bSJed Brown   PetscScalar       *f;
546c4762a1bSJed Brown   const PetscScalar *y;
547c4762a1bSJed Brown 
5487510d9b0SBarry Smith   PetscFunctionBeginUser;
5499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
5509566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
551c4762a1bSJed Brown   f[0] = -y[1] - y[0] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
552c4762a1bSJed Brown   f[1] = y[0] - y[1] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
553c4762a1bSJed Brown   f[2] = y[0] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
5549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
5559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
556c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
5579566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
558c4762a1bSJed Brown   PetscFunctionReturn(0);
559c4762a1bSJed Brown }
560c4762a1bSJed Brown 
561*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
562c4762a1bSJed Brown   const PetscScalar *y;
563c4762a1bSJed Brown   PetscInt           row[3] = {0, 1, 2};
564c4762a1bSJed Brown   PetscScalar        value[3][3], fac, fac2;
565c4762a1bSJed Brown 
5667510d9b0SBarry Smith   PetscFunctionBeginUser;
5679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
568c4762a1bSJed Brown   fac         = PetscPowScalar(y[0] * y[0] + y[1] * y[1], -1.5);
569c4762a1bSJed Brown   fac2        = PetscPowScalar(y[0] * y[0] + y[1] * y[1], -0.5);
570c4762a1bSJed Brown   value[0][0] = a + (y[1] * y[1] * y[2]) * fac;
571c4762a1bSJed Brown   value[0][1] = 1.0 - (y[0] * y[1] * y[2]) * fac;
572c4762a1bSJed Brown   value[0][2] = y[0] * fac2;
573c4762a1bSJed Brown   value[1][0] = -1.0 - y[0] * y[1] * y[2] * fac;
574c4762a1bSJed Brown   value[1][1] = a + y[0] * y[0] * y[2] * fac;
575c4762a1bSJed Brown   value[1][2] = y[1] * fac2;
576c4762a1bSJed Brown   value[2][0] = -y[1] * y[1] * fac;
577c4762a1bSJed Brown   value[2][1] = y[0] * y[1] * fac;
578c4762a1bSJed Brown   value[2][2] = a;
5799566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES));
5809566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
583c4762a1bSJed Brown   PetscFunctionReturn(0);
584c4762a1bSJed Brown }
585c4762a1bSJed Brown 
586c4762a1bSJed Brown /* Hull, 1972, Problem B5 */
587c4762a1bSJed Brown 
588*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
589c4762a1bSJed Brown   PetscScalar       *f;
590c4762a1bSJed Brown   const PetscScalar *y;
591c4762a1bSJed Brown 
5927510d9b0SBarry Smith   PetscFunctionBeginUser;
5939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
5949566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
595c4762a1bSJed Brown   f[0] = y[1] * y[2];
596c4762a1bSJed Brown   f[1] = -y[0] * y[2];
597c4762a1bSJed Brown   f[2] = -0.51 * y[0] * y[1];
5989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
5999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
600c4762a1bSJed Brown   PetscFunctionReturn(0);
601c4762a1bSJed Brown }
602c4762a1bSJed Brown 
603*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
604c4762a1bSJed Brown   PetscScalar       *f;
605c4762a1bSJed Brown   const PetscScalar *y;
606c4762a1bSJed Brown 
6077510d9b0SBarry Smith   PetscFunctionBeginUser;
6089566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
6099566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
610c4762a1bSJed Brown   f[0] = y[1] * y[2];
611c4762a1bSJed Brown   f[1] = -y[0] * y[2];
612c4762a1bSJed Brown   f[2] = -0.51 * y[0] * y[1];
6139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
6149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
615c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
6169566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
617c4762a1bSJed Brown   PetscFunctionReturn(0);
618c4762a1bSJed Brown }
619c4762a1bSJed Brown 
620*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
621c4762a1bSJed Brown   const PetscScalar *y;
622c4762a1bSJed Brown   PetscInt           row[3] = {0, 1, 2};
623c4762a1bSJed Brown   PetscScalar        value[3][3];
624c4762a1bSJed Brown 
6257510d9b0SBarry Smith   PetscFunctionBeginUser;
6269566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
627*9371c9d4SSatish Balay   value[0][0] = a;
628*9371c9d4SSatish Balay   value[0][1] = -y[2];
629*9371c9d4SSatish Balay   value[0][2] = -y[1];
630*9371c9d4SSatish Balay   value[1][0] = y[2];
631*9371c9d4SSatish Balay   value[1][1] = a;
632*9371c9d4SSatish Balay   value[1][2] = y[0];
633*9371c9d4SSatish Balay   value[2][0] = 0.51 * y[1];
634*9371c9d4SSatish Balay   value[2][1] = 0.51 * y[0];
635*9371c9d4SSatish Balay   value[2][2] = a;
6369566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES));
6379566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
6389566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
6399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
640c4762a1bSJed Brown   PetscFunctionReturn(0);
641c4762a1bSJed Brown }
642c4762a1bSJed Brown 
643c4762a1bSJed Brown /* Kulikov, 2013, Problem I */
644c4762a1bSJed Brown 
645*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
646c4762a1bSJed Brown   PetscScalar       *f;
647c4762a1bSJed Brown   const PetscScalar *y;
648c4762a1bSJed Brown 
6497510d9b0SBarry Smith   PetscFunctionBeginUser;
6509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
6519566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
652c4762a1bSJed Brown   f[0] = 2. * t * PetscPowScalar(y[1], 1. / 5.) * y[3];
653c4762a1bSJed Brown   f[1] = 10. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.));
654c4762a1bSJed Brown   f[2] = 2. * t * y[3];
655c4762a1bSJed Brown   f[3] = -2. * t * PetscLogScalar(y[0]);
6569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
6579566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
658c4762a1bSJed Brown   PetscFunctionReturn(0);
659c4762a1bSJed Brown }
660c4762a1bSJed Brown 
661*9371c9d4SSatish Balay PetscErrorCode RHSJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s) {
662c4762a1bSJed Brown   const PetscScalar *y;
663c4762a1bSJed Brown   PetscInt           row[4] = {0, 1, 2, 3};
664c4762a1bSJed Brown   PetscScalar        value[4][4];
665c4762a1bSJed Brown   PetscScalar        m1, m2;
6667510d9b0SBarry Smith   PetscFunctionBeginUser;
6679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
668c4762a1bSJed Brown   m1          = (2. * t * y[3]) / (5. * PetscPowScalar(y[1], 4. / 5.));
669c4762a1bSJed Brown   m2          = 2. * t * PetscPowScalar(y[1], 1. / 5.);
670*9371c9d4SSatish Balay   value[0][0] = 0.;
671*9371c9d4SSatish Balay   value[0][1] = m1;
672*9371c9d4SSatish Balay   value[0][2] = 0.;
673*9371c9d4SSatish Balay   value[0][3] = m2;
674c4762a1bSJed Brown   m1          = 50. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.));
675c4762a1bSJed Brown   m2          = 10. * t * PetscExpScalar(5.0 * (y[2] - 1.));
676*9371c9d4SSatish Balay   value[1][0] = 0.;
677*9371c9d4SSatish Balay   value[1][1] = 0.;
678*9371c9d4SSatish Balay   value[1][2] = m1;
679*9371c9d4SSatish Balay   value[1][3] = m2;
680*9371c9d4SSatish Balay   value[2][0] = 0.;
681*9371c9d4SSatish Balay   value[2][1] = 0.;
682*9371c9d4SSatish Balay   value[2][2] = 0.;
683*9371c9d4SSatish Balay   value[2][3] = 2 * t;
684*9371c9d4SSatish Balay   value[3][0] = -2. * t / y[0];
685*9371c9d4SSatish Balay   value[3][1] = 0.;
686*9371c9d4SSatish Balay   value[3][2] = 0.;
687*9371c9d4SSatish Balay   value[3][3] = 0.;
6889566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 4, &row[0], 4, &row[0], &value[0][0], INSERT_VALUES));
6899566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
6909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
6919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
692c4762a1bSJed Brown   PetscFunctionReturn(0);
693c4762a1bSJed Brown }
694c4762a1bSJed Brown 
695*9371c9d4SSatish Balay PetscErrorCode IFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
696c4762a1bSJed Brown   PetscScalar       *f;
697c4762a1bSJed Brown   const PetscScalar *y;
698c4762a1bSJed Brown 
6997510d9b0SBarry Smith   PetscFunctionBeginUser;
7009566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
7019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
702c4762a1bSJed Brown   f[0] = 2. * t * PetscPowScalar(y[1], 1. / 5.) * y[3];
703c4762a1bSJed Brown   f[1] = 10. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.));
704c4762a1bSJed Brown   f[2] = 2. * t * y[3];
705c4762a1bSJed Brown   f[3] = -2. * t * PetscLogScalar(y[0]);
7069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
7079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
708c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
7099566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
710c4762a1bSJed Brown   PetscFunctionReturn(0);
711c4762a1bSJed Brown }
712c4762a1bSJed Brown 
713*9371c9d4SSatish Balay PetscErrorCode IJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
714c4762a1bSJed Brown   const PetscScalar *y;
715c4762a1bSJed Brown   PetscInt           row[4] = {0, 1, 2, 3};
716c4762a1bSJed Brown   PetscScalar        value[4][4];
717c4762a1bSJed Brown   PetscScalar        m1, m2;
718c4762a1bSJed Brown 
7197510d9b0SBarry Smith   PetscFunctionBeginUser;
7209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
721c4762a1bSJed Brown   m1          = (2. * t * y[3]) / (5. * PetscPowScalar(y[1], 4. / 5.));
722c4762a1bSJed Brown   m2          = 2. * t * PetscPowScalar(y[1], 1. / 5.);
723*9371c9d4SSatish Balay   value[0][0] = a;
724*9371c9d4SSatish Balay   value[0][1] = m1;
725*9371c9d4SSatish Balay   value[0][2] = 0.;
726*9371c9d4SSatish Balay   value[0][3] = m2;
727c4762a1bSJed Brown   m1          = 50. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.));
728c4762a1bSJed Brown   m2          = 10. * t * PetscExpScalar(5.0 * (y[2] - 1.));
729*9371c9d4SSatish Balay   value[1][0] = 0.;
730*9371c9d4SSatish Balay   value[1][1] = a;
731*9371c9d4SSatish Balay   value[1][2] = m1;
732*9371c9d4SSatish Balay   value[1][3] = m2;
733*9371c9d4SSatish Balay   value[2][0] = 0.;
734*9371c9d4SSatish Balay   value[2][1] = 0.;
735*9371c9d4SSatish Balay   value[2][2] = a;
736*9371c9d4SSatish Balay   value[2][3] = 2 * t;
737*9371c9d4SSatish Balay   value[3][0] = -2. * t / y[0];
738*9371c9d4SSatish Balay   value[3][1] = 0.;
739*9371c9d4SSatish Balay   value[3][2] = 0.;
740*9371c9d4SSatish Balay   value[3][3] = a;
7419566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 4, &row[0], 4, &row[0], &value[0][0], INSERT_VALUES));
7429566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7439566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
7449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
745c4762a1bSJed Brown   PetscFunctionReturn(0);
746c4762a1bSJed Brown }
747c4762a1bSJed Brown 
748c4762a1bSJed Brown /* Hull, 1972, Problem C1 */
749c4762a1bSJed Brown 
750*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
751c4762a1bSJed Brown   PetscScalar       *f;
752c4762a1bSJed Brown   const PetscScalar *y;
753c4762a1bSJed Brown   PetscInt           N, i;
754c4762a1bSJed Brown 
7557510d9b0SBarry Smith   PetscFunctionBeginUser;
7569566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
7579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
7589566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
759c4762a1bSJed Brown   f[0] = -y[0];
760*9371c9d4SSatish Balay   for (i = 1; i < N - 1; i++) { f[i] = y[i - 1] - y[i]; }
761c4762a1bSJed Brown   f[N - 1] = y[N - 2];
7629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
7639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
764c4762a1bSJed Brown   PetscFunctionReturn(0);
765c4762a1bSJed Brown }
766c4762a1bSJed Brown 
767*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
768c4762a1bSJed Brown   PetscScalar       *f;
769c4762a1bSJed Brown   const PetscScalar *y;
770c4762a1bSJed Brown   PetscInt           N, i;
771c4762a1bSJed Brown 
7727510d9b0SBarry Smith   PetscFunctionBeginUser;
7739566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
7749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
7759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
776c4762a1bSJed Brown   f[0] = -y[0];
777*9371c9d4SSatish Balay   for (i = 1; i < N - 1; i++) { f[i] = y[i - 1] - y[i]; }
778c4762a1bSJed Brown   f[N - 1] = y[N - 2];
7799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
7809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
781c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
7829566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
783c4762a1bSJed Brown   PetscFunctionReturn(0);
784c4762a1bSJed Brown }
785c4762a1bSJed Brown 
786*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
787c4762a1bSJed Brown   const PetscScalar *y;
788c4762a1bSJed Brown   PetscInt           N, i, col[2];
789c4762a1bSJed Brown   PetscScalar        value[2];
790c4762a1bSJed Brown 
7917510d9b0SBarry Smith   PetscFunctionBeginUser;
7929566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
7939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
794c4762a1bSJed Brown   i        = 0;
795*9371c9d4SSatish Balay   value[0] = a + 1;
796*9371c9d4SSatish Balay   col[0]   = 0;
797*9371c9d4SSatish Balay   value[1] = 0;
798*9371c9d4SSatish Balay   col[1]   = 1;
7999566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
800c4762a1bSJed Brown   for (i = 0; i < N; i++) {
801*9371c9d4SSatish Balay     value[0] = -1;
802*9371c9d4SSatish Balay     col[0]   = i - 1;
803*9371c9d4SSatish Balay     value[1] = a + 1;
804*9371c9d4SSatish Balay     col[1]   = i;
8059566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
806c4762a1bSJed Brown   }
807c4762a1bSJed Brown   i        = N - 1;
808*9371c9d4SSatish Balay   value[0] = -1;
809*9371c9d4SSatish Balay   col[0]   = N - 2;
810*9371c9d4SSatish Balay   value[1] = a;
811*9371c9d4SSatish Balay   col[1]   = N - 1;
8129566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
8139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
8149566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
8159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
816c4762a1bSJed Brown   PetscFunctionReturn(0);
817c4762a1bSJed Brown }
818c4762a1bSJed Brown 
819c4762a1bSJed Brown /* Hull, 1972, Problem C2 */
820c4762a1bSJed Brown 
821*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
822c4762a1bSJed Brown   const PetscScalar *y;
823c4762a1bSJed Brown   PetscScalar       *f;
824c4762a1bSJed Brown   PetscInt           N, i;
825c4762a1bSJed Brown 
8267510d9b0SBarry Smith   PetscFunctionBeginUser;
8279566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
8299566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
830c4762a1bSJed Brown   f[0] = -y[0];
831*9371c9d4SSatish Balay   for (i = 1; i < N - 1; i++) { f[i] = (PetscReal)i * y[i - 1] - (PetscReal)(i + 1) * y[i]; }
832c4762a1bSJed Brown   f[N - 1] = (PetscReal)(N - 1) * y[N - 2];
8339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
8349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
835c4762a1bSJed Brown   PetscFunctionReturn(0);
836c4762a1bSJed Brown }
837c4762a1bSJed Brown 
838*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
839c4762a1bSJed Brown   PetscScalar       *f;
840c4762a1bSJed Brown   const PetscScalar *y;
841c4762a1bSJed Brown   PetscInt           N, i;
842c4762a1bSJed Brown 
8437510d9b0SBarry Smith   PetscFunctionBeginUser;
8449566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
8469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
847c4762a1bSJed Brown   f[0] = -y[0];
848*9371c9d4SSatish Balay   for (i = 1; i < N - 1; i++) { f[i] = (PetscReal)i * y[i - 1] - (PetscReal)(i + 1) * y[i]; }
849c4762a1bSJed Brown   f[N - 1] = (PetscReal)(N - 1) * y[N - 2];
8509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
8519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
852c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
8539566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
854c4762a1bSJed Brown   PetscFunctionReturn(0);
855c4762a1bSJed Brown }
856c4762a1bSJed Brown 
857*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
858c4762a1bSJed Brown   const PetscScalar *y;
859c4762a1bSJed Brown   PetscInt           N, i, col[2];
860c4762a1bSJed Brown   PetscScalar        value[2];
861c4762a1bSJed Brown 
8627510d9b0SBarry Smith   PetscFunctionBeginUser;
8639566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
865c4762a1bSJed Brown   i        = 0;
866*9371c9d4SSatish Balay   value[0] = a + 1;
867*9371c9d4SSatish Balay   col[0]   = 0;
868*9371c9d4SSatish Balay   value[1] = 0;
869*9371c9d4SSatish Balay   col[1]   = 1;
8709566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
871c4762a1bSJed Brown   for (i = 0; i < N; i++) {
872*9371c9d4SSatish Balay     value[0] = -(PetscReal)i;
873*9371c9d4SSatish Balay     col[0]   = i - 1;
874*9371c9d4SSatish Balay     value[1] = a + (PetscReal)(i + 1);
875*9371c9d4SSatish Balay     col[1]   = i;
8769566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
877c4762a1bSJed Brown   }
878c4762a1bSJed Brown   i        = N - 1;
879*9371c9d4SSatish Balay   value[0] = -(PetscReal)(N - 1);
880*9371c9d4SSatish Balay   col[0]   = N - 2;
881*9371c9d4SSatish Balay   value[1] = a;
882*9371c9d4SSatish Balay   col[1]   = N - 1;
8839566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
8849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
8859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
8869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
887c4762a1bSJed Brown   PetscFunctionReturn(0);
888c4762a1bSJed Brown }
889c4762a1bSJed Brown 
890c4762a1bSJed Brown /* Hull, 1972, Problem C3 and C4 */
891c4762a1bSJed Brown 
892*9371c9d4SSatish Balay PetscErrorCode RHSFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec F, void *s) {
893c4762a1bSJed Brown   PetscScalar       *f;
894c4762a1bSJed Brown   const PetscScalar *y;
895c4762a1bSJed Brown   PetscInt           N, i;
896c4762a1bSJed Brown 
8977510d9b0SBarry Smith   PetscFunctionBeginUser;
8989566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
9009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
901c4762a1bSJed Brown   f[0] = -2.0 * y[0] + y[1];
902*9371c9d4SSatish Balay   for (i = 1; i < N - 1; i++) { f[i] = y[i - 1] - 2.0 * y[i] + y[i + 1]; }
903c4762a1bSJed Brown   f[N - 1] = y[N - 2] - 2.0 * y[N - 1];
9049566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
9059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
906c4762a1bSJed Brown   PetscFunctionReturn(0);
907c4762a1bSJed Brown }
908c4762a1bSJed Brown 
909*9371c9d4SSatish Balay PetscErrorCode IFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s) {
910c4762a1bSJed Brown   PetscScalar       *f;
911c4762a1bSJed Brown   const PetscScalar *y;
912c4762a1bSJed Brown   PetscInt           N, i;
913c4762a1bSJed Brown 
9147510d9b0SBarry Smith   PetscFunctionBeginUser;
9159566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
9169566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
9179566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
918c4762a1bSJed Brown   f[0] = -2.0 * y[0] + y[1];
919*9371c9d4SSatish Balay   for (i = 1; i < N - 1; i++) { f[i] = y[i - 1] - 2.0 * y[i] + y[i + 1]; }
920c4762a1bSJed Brown   f[N - 1] = y[N - 2] - 2.0 * y[N - 1];
9219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
9229566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
923c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
9249566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
925c4762a1bSJed Brown   PetscFunctionReturn(0);
926c4762a1bSJed Brown }
927c4762a1bSJed Brown 
928*9371c9d4SSatish Balay PetscErrorCode IJacobian_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s) {
929c4762a1bSJed Brown   const PetscScalar *y;
930c4762a1bSJed Brown   PetscScalar        value[3];
931c4762a1bSJed Brown   PetscInt           N, i, col[3];
932c4762a1bSJed Brown 
9337510d9b0SBarry Smith   PetscFunctionBeginUser;
9349566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
9359566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
936c4762a1bSJed Brown   for (i = 0; i < N; i++) {
937c4762a1bSJed Brown     if (i == 0) {
938*9371c9d4SSatish Balay       value[0] = a + 2;
939*9371c9d4SSatish Balay       col[0]   = i;
940*9371c9d4SSatish Balay       value[1] = -1;
941*9371c9d4SSatish Balay       col[1]   = i + 1;
942*9371c9d4SSatish Balay       value[2] = 0;
943*9371c9d4SSatish Balay       col[2]   = i + 2;
944c4762a1bSJed Brown     } else if (i == N - 1) {
945*9371c9d4SSatish Balay       value[0] = 0;
946*9371c9d4SSatish Balay       col[0]   = i - 2;
947*9371c9d4SSatish Balay       value[1] = -1;
948*9371c9d4SSatish Balay       col[1]   = i - 1;
949*9371c9d4SSatish Balay       value[2] = a + 2;
950*9371c9d4SSatish Balay       col[2]   = i;
951c4762a1bSJed Brown     } else {
952*9371c9d4SSatish Balay       value[0] = -1;
953*9371c9d4SSatish Balay       col[0]   = i - 1;
954*9371c9d4SSatish Balay       value[1] = a + 2;
955*9371c9d4SSatish Balay       col[1]   = i;
956*9371c9d4SSatish Balay       value[2] = -1;
957*9371c9d4SSatish Balay       col[2]   = i + 1;
958c4762a1bSJed Brown     }
9599566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
960c4762a1bSJed Brown   }
9619566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
9629566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
9639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
964c4762a1bSJed Brown   PetscFunctionReturn(0);
965c4762a1bSJed Brown }
966c4762a1bSJed Brown 
967c4762a1bSJed Brown /***************************************************************************/
968c4762a1bSJed Brown 
969c4762a1bSJed Brown /* Sets the initial solution for the IVP and sets up the function pointers*/
970*9371c9d4SSatish Balay PetscErrorCode Initialize(Vec Y, void *s) {
971c4762a1bSJed Brown   char        *p = (char *)s;
972c4762a1bSJed Brown   PetscScalar *y;
973c4762a1bSJed Brown   PetscReal    t0;
974c4762a1bSJed Brown   PetscInt     N = GetSize((const char *)s);
975c4762a1bSJed Brown   PetscBool    flg;
976c4762a1bSJed Brown 
9777510d9b0SBarry Smith   PetscFunctionBeginUser;
978c4762a1bSJed Brown   VecZeroEntries(Y);
9799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y, &y));
980c4762a1bSJed Brown   if (!strcmp(p, "hull1972a1")) {
981c4762a1bSJed Brown     y[0]        = 1.0;
982c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A1;
983c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A1;
984c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A1;
985c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A1;
986c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a2")) {
987c4762a1bSJed Brown     y[0]        = 1.0;
988c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A2;
989c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A2;
990c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A2;
991c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A2;
992c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a3")) {
993c4762a1bSJed Brown     y[0]        = 1.0;
994c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A3;
995c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A3;
996c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A3;
997c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A3;
998c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a4")) {
999c4762a1bSJed Brown     y[0]        = 1.0;
1000c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A4;
1001c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A4;
1002c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A4;
1003c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A4;
1004c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a5")) {
1005c4762a1bSJed Brown     y[0]        = 4.0;
1006c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A5;
1007c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A5;
1008c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A5;
1009c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A5;
1010c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b1")) {
1011c4762a1bSJed Brown     y[0]        = 1.0;
1012c4762a1bSJed Brown     y[1]        = 3.0;
1013c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B1;
1014c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B1;
1015c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B1;
1016c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b2")) {
1017c4762a1bSJed Brown     y[0]        = 2.0;
1018c4762a1bSJed Brown     y[1]        = 0.0;
1019c4762a1bSJed Brown     y[2]        = 1.0;
1020c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B2;
1021c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B2;
1022c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B2;
1023c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b3")) {
1024c4762a1bSJed Brown     y[0]        = 1.0;
1025c4762a1bSJed Brown     y[1]        = 0.0;
1026c4762a1bSJed Brown     y[2]        = 0.0;
1027c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B3;
1028c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B3;
1029c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B3;
1030c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b4")) {
1031c4762a1bSJed Brown     y[0]        = 3.0;
1032c4762a1bSJed Brown     y[1]        = 0.0;
1033c4762a1bSJed Brown     y[2]        = 0.0;
1034c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B4;
1035c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B4;
1036c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B4;
1037c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b5")) {
1038c4762a1bSJed Brown     y[0]        = 0.0;
1039c4762a1bSJed Brown     y[1]        = 1.0;
1040c4762a1bSJed Brown     y[2]        = 1.0;
1041c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B5;
1042c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B5;
1043c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B5;
1044c4762a1bSJed Brown   } else if (!strcmp(p, "kulik2013i")) {
1045c4762a1bSJed Brown     t0          = 0.;
1046c4762a1bSJed Brown     y[0]        = PetscExpReal(PetscSinReal(t0 * t0));
1047c4762a1bSJed Brown     y[1]        = PetscExpReal(5. * PetscSinReal(t0 * t0));
1048c4762a1bSJed Brown     y[2]        = PetscSinReal(t0 * t0) + 1.0;
1049c4762a1bSJed Brown     y[3]        = PetscCosReal(t0 * t0);
1050c4762a1bSJed Brown     RHSFunction = RHSFunction_Kulikov2013I;
1051c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Kulikov2013I;
1052c4762a1bSJed Brown     IFunction   = IFunction_Kulikov2013I;
1053c4762a1bSJed Brown     IJacobian   = IJacobian_Kulikov2013I;
1054c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972c1")) {
1055c4762a1bSJed Brown     y[0]        = 1.0;
1056c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C1;
1057c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C1;
1058c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C1;
1059c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972c2")) {
1060c4762a1bSJed Brown     y[0]        = 1.0;
1061c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C2;
1062c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C2;
1063c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C2;
1064*9371c9d4SSatish Balay   } else if ((!strcmp(p, "hull1972c3")) || (!strcmp(p, "hull1972c4"))) {
1065c4762a1bSJed Brown     y[0]        = 1.0;
1066c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C34;
1067c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C34;
1068c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C34;
1069c4762a1bSJed Brown   }
10709566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalarArray(NULL, NULL, "-yinit", y, &N, &flg));
107163a3b9bcSJacob Faibussowitsch   PetscCheck((N == GetSize((const char *)s)) || !flg, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Number of initial values %" PetscInt_FMT " does not match problem size %" PetscInt_FMT ".", N, GetSize((const char *)s));
10729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y, &y));
1073c4762a1bSJed Brown   PetscFunctionReturn(0);
1074c4762a1bSJed Brown }
1075c4762a1bSJed Brown 
1076c4762a1bSJed Brown /* Calculates the exact solution to problems that have one */
1077*9371c9d4SSatish Balay PetscErrorCode ExactSolution(Vec Y, void *s, PetscReal t, PetscBool *flag) {
1078c4762a1bSJed Brown   char        *p = (char *)s;
1079c4762a1bSJed Brown   PetscScalar *y;
1080c4762a1bSJed Brown 
10817510d9b0SBarry Smith   PetscFunctionBeginUser;
1082c4762a1bSJed Brown   if (!strcmp(p, "hull1972a1")) {
10839566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1084c4762a1bSJed Brown     y[0]  = PetscExpReal(-t);
1085c4762a1bSJed Brown     *flag = PETSC_TRUE;
10869566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1087c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a2")) {
10889566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1089c4762a1bSJed Brown     y[0]  = 1.0 / PetscSqrtReal(t + 1);
1090c4762a1bSJed Brown     *flag = PETSC_TRUE;
10919566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1092c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a3")) {
10939566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1094c4762a1bSJed Brown     y[0]  = PetscExpReal(PetscSinReal(t));
1095c4762a1bSJed Brown     *flag = PETSC_TRUE;
10969566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1097c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a4")) {
10989566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1099c4762a1bSJed Brown     y[0]  = 20.0 / (1 + 19.0 * PetscExpReal(-t / 4.0));
1100c4762a1bSJed Brown     *flag = PETSC_TRUE;
11019566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1102c4762a1bSJed Brown   } else if (!strcmp(p, "kulik2013i")) {
11039566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1104c4762a1bSJed Brown     y[0]  = PetscExpReal(PetscSinReal(t * t));
1105c4762a1bSJed Brown     y[1]  = PetscExpReal(5. * PetscSinReal(t * t));
1106c4762a1bSJed Brown     y[2]  = PetscSinReal(t * t) + 1.0;
1107c4762a1bSJed Brown     y[3]  = PetscCosReal(t * t);
1108c4762a1bSJed Brown     *flag = PETSC_TRUE;
11099566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1110c4762a1bSJed Brown   } else {
11119566063dSJacob Faibussowitsch     PetscCall(VecSet(Y, 0));
1112c4762a1bSJed Brown     *flag = PETSC_FALSE;
1113c4762a1bSJed Brown   }
1114c4762a1bSJed Brown   PetscFunctionReturn(0);
1115c4762a1bSJed Brown }
1116c4762a1bSJed Brown 
1117c4762a1bSJed Brown /* Solves the specified ODE and computes the error if exact solution is available */
1118*9371c9d4SSatish Balay PetscErrorCode SolveODE(char *ptype, PetscReal dt, PetscReal tfinal, PetscInt maxiter, PetscReal *error, PetscBool *exact_flag) {
1119c4762a1bSJed Brown   TS        ts;          /* time-integrator                        */
1120c4762a1bSJed Brown   Vec       Y;           /* Solution vector                        */
1121c4762a1bSJed Brown   Vec       Yex;         /* Exact solution                         */
1122c4762a1bSJed Brown   PetscInt  N;           /* Size of the system of equations        */
1123c4762a1bSJed Brown   TSType    time_scheme; /* Type of time-integration scheme        */
1124c4762a1bSJed Brown   Mat       Jac = NULL;  /* Jacobian matrix                        */
1125c4762a1bSJed Brown   Vec       Yerr;        /* Auxiliary solution vector              */
1126c4762a1bSJed Brown   PetscReal err_norm;    /* Estimated error norm                   */
1127c4762a1bSJed Brown   PetscReal final_time;  /* Actual final time from the integrator  */
1128c4762a1bSJed Brown 
11297510d9b0SBarry Smith   PetscFunctionBeginUser;
1130c4762a1bSJed Brown   N = GetSize((const char *)&ptype[0]);
11313c633725SBarry Smith   PetscCheck(N >= 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Illegal problem specification.");
11329566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &Y));
11339566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(Y, N, PETSC_DECIDE));
11349566063dSJacob Faibussowitsch   PetscCall(VecSetUp(Y));
11359566063dSJacob Faibussowitsch   PetscCall(VecSet(Y, 0));
1136c4762a1bSJed Brown 
1137c4762a1bSJed Brown   /* Initialize the problem */
11389566063dSJacob Faibussowitsch   PetscCall(Initialize(Y, &ptype[0]));
1139c4762a1bSJed Brown 
1140c4762a1bSJed Brown   /* Create and initialize the time-integrator                            */
11419566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1142c4762a1bSJed Brown   /* Default time integration options                                     */
11439566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSRK));
11449566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, maxiter));
11459566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, tfinal));
11469566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
11479566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1148c4762a1bSJed Brown   /* Read command line options for time integration                       */
11499566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1150c4762a1bSJed Brown   /* Set solution vector                                                  */
11519566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, Y));
1152c4762a1bSJed Brown   /* Specify left/right-hand side functions                               */
11539566063dSJacob Faibussowitsch   PetscCall(TSGetType(ts, &time_scheme));
1154c4762a1bSJed Brown 
1155c4762a1bSJed Brown   if ((!strcmp(time_scheme, TSEULER)) || (!strcmp(time_scheme, TSRK)) || (!strcmp(time_scheme, TSSSP) || (!strcmp(time_scheme, TSGLEE)))) {
1156c4762a1bSJed Brown     /* Explicit time-integration -> specify right-hand side function ydot = f(y) */
11579566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &ptype[0]));
11589566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &Jac));
11599566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Jac, PETSC_DECIDE, PETSC_DECIDE, N, N));
11609566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(Jac));
11619566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Jac));
11629566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, Jac, Jac, RHSJacobian, &ptype[0]));
1163c4762a1bSJed Brown   } else if ((!strcmp(time_scheme, TSTHETA)) || (!strcmp(time_scheme, TSBEULER)) || (!strcmp(time_scheme, TSCN)) || (!strcmp(time_scheme, TSALPHA)) || (!strcmp(time_scheme, TSARKIMEX))) {
1164c4762a1bSJed Brown     /* Implicit time-integration -> specify left-hand side function ydot-f(y) = 0 */
1165c4762a1bSJed Brown     /* and its Jacobian function                                                 */
11669566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction, &ptype[0]));
11679566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &Jac));
11689566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Jac, PETSC_DECIDE, PETSC_DECIDE, N, N));
11699566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(Jac));
11709566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Jac));
11719566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, Jac, Jac, IJacobian, &ptype[0]));
1172c4762a1bSJed Brown   }
1173c4762a1bSJed Brown 
1174c4762a1bSJed Brown   /* Solve */
11759566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, Y));
11769566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &final_time));
1177c4762a1bSJed Brown 
1178c4762a1bSJed Brown   /* Get the estimated error, if available */
11799566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(Y, &Yerr));
11809566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Yerr));
11819566063dSJacob Faibussowitsch   PetscCall(TSGetTimeError(ts, 0, &Yerr));
11829566063dSJacob Faibussowitsch   PetscCall(VecNorm(Yerr, NORM_2, &err_norm));
11839566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Yerr));
118463a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Estimated Error = %e.\n", (double)err_norm));
1185c4762a1bSJed Brown 
1186c4762a1bSJed Brown   /* Exact solution */
11879566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(Y, &Yex));
118863a3b9bcSJacob Faibussowitsch   if (PetscAbsReal(final_time - tfinal) > 2. * PETSC_MACHINE_EPSILON) {
11899566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Note: There is a difference between the prescribed final time %g and the actual final time, %g.\n", (double)tfinal, (double)final_time));
1190c4762a1bSJed Brown   }
11919566063dSJacob Faibussowitsch   PetscCall(ExactSolution(Yex, &ptype[0], final_time, exact_flag));
1192c4762a1bSJed Brown 
1193c4762a1bSJed Brown   /* Calculate Error */
11949566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Yex, -1.0, Y));
11959566063dSJacob Faibussowitsch   PetscCall(VecNorm(Yex, NORM_2, error));
1196c4762a1bSJed Brown   *error = PetscSqrtReal(((*error) * (*error)) / N);
1197c4762a1bSJed Brown 
1198c4762a1bSJed Brown   /* Clean up and finalize */
11999566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jac));
12009566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
12019566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Yex));
12029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
1203c4762a1bSJed Brown 
1204c4762a1bSJed Brown   PetscFunctionReturn(0);
1205c4762a1bSJed Brown }
1206c4762a1bSJed Brown 
1207*9371c9d4SSatish Balay int main(int argc, char **argv) {
1208c4762a1bSJed Brown   char        ptype[256] = "hull1972a1"; /* Problem specification                                */
1209c4762a1bSJed Brown   PetscInt    n_refine   = 1;            /* Number of refinement levels for convergence analysis */
1210c4762a1bSJed Brown   PetscReal   refine_fac = 2.0;          /* Refinement factor for dt                             */
1211c4762a1bSJed Brown   PetscReal   dt_initial = 0.01;         /* Initial default value of dt                          */
1212c4762a1bSJed Brown   PetscReal   dt;
1213c4762a1bSJed Brown   PetscReal   tfinal  = 20.0;   /* Final time for the time-integration                  */
1214c4762a1bSJed Brown   PetscInt    maxiter = 100000; /* Maximum number of time-integration iterations        */
1215c4762a1bSJed Brown   PetscReal  *error;            /* Array to store the errors for convergence analysis   */
1216c4762a1bSJed Brown   PetscMPIInt size;             /* No of processors                                     */
1217c4762a1bSJed Brown   PetscBool   flag;             /* Flag denoting availability of exact solution         */
1218c4762a1bSJed Brown   PetscInt    r;
1219c4762a1bSJed Brown 
1220c4762a1bSJed Brown   /* Initialize program */
1221327415f7SBarry Smith   PetscFunctionBeginUser;
12229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1223c4762a1bSJed Brown 
1224c4762a1bSJed Brown   /* Check if running with only 1 proc */
12259566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
12263c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
1227c4762a1bSJed Brown 
1228d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex31", NULL);
12299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-problem", "Problem specification", "<hull1972a1>", ptype, ptype, sizeof(ptype), NULL));
12309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-refinement_levels", "Number of refinement levels for convergence analysis", "<1>", n_refine, &n_refine, NULL));
12319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-refinement_factor", "Refinement factor for dt", "<2.0>", refine_fac, &refine_fac, NULL));
12329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dt", "Time step size (for convergence analysis, initial time step)", "<0.01>", dt_initial, &dt_initial, NULL));
12339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-final_time", "Final time for the time-integration", "<20.0>", tfinal, &tfinal, NULL));
1234d0609cedSBarry Smith   PetscOptionsEnd();
1235c4762a1bSJed Brown 
12369566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_refine, &error));
1237c4762a1bSJed Brown   for (r = 0, dt = dt_initial; r < n_refine; r++) {
1238c4762a1bSJed Brown     error[r] = 0;
1239c4762a1bSJed Brown     if (r > 0) dt /= refine_fac;
1240c4762a1bSJed Brown 
124163a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solving ODE \"%s\" with dt %f, final time %f and system size %" PetscInt_FMT ".\n", ptype, (double)dt, (double)tfinal, GetSize(&ptype[0])));
12429566063dSJacob Faibussowitsch     PetscCall(SolveODE(&ptype[0], dt, tfinal, maxiter, &error[r], &flag));
1243c4762a1bSJed Brown     if (flag) {
1244c4762a1bSJed Brown       /* If exact solution available for the specified ODE */
1245c4762a1bSJed Brown       if (r > 0) {
1246c4762a1bSJed Brown         PetscReal conv_rate = (PetscLogReal(error[r]) - PetscLogReal(error[r - 1])) / (-PetscLogReal(refine_fac));
12479566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error           = %E,\tConvergence rate = %f.\n", (double)error[r], (double)conv_rate));
1248c4762a1bSJed Brown       } else {
124963a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error           = %E.\n", (double)error[r]));
1250c4762a1bSJed Brown       }
1251c4762a1bSJed Brown     }
1252c4762a1bSJed Brown   }
12539566063dSJacob Faibussowitsch   PetscCall(PetscFree(error));
12549566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1255b122ec5aSJacob Faibussowitsch   return 0;
1256c4762a1bSJed Brown }
1257c4762a1bSJed Brown 
1258c4762a1bSJed Brown /*TEST
1259c4762a1bSJed Brown 
1260c4762a1bSJed Brown     test:
1261c4762a1bSJed Brown       suffix: 2
1262c4762a1bSJed Brown       args: -ts_type glee -final_time 5 -ts_adapt_type none
1263c4762a1bSJed Brown       timeoutfactor: 3
1264c4762a1bSJed Brown       requires: !single
1265c4762a1bSJed Brown 
1266c4762a1bSJed Brown     test:
1267c4762a1bSJed Brown       suffix: 3
1268c4762a1bSJed Brown       args: -ts_type glee -final_time 5 -ts_adapt_type glee -ts_adapt_monitor  -ts_max_steps 50  -problem hull1972a3 -ts_adapt_glee_use_local 1
1269c4762a1bSJed Brown       timeoutfactor: 3
1270c4762a1bSJed Brown       requires: !single
1271c4762a1bSJed Brown 
1272c4762a1bSJed Brown     test:
1273c4762a1bSJed Brown       suffix: 4
1274c4762a1bSJed Brown       args: -ts_type glee -final_time 5 -ts_adapt_type glee -ts_adapt_monitor  -ts_max_steps 50  -problem hull1972a3  -ts_max_reject 100 -ts_adapt_glee_use_local 0
1275c4762a1bSJed Brown       timeoutfactor: 3
1276c4762a1bSJed Brown       requires: !single !__float128
1277c4762a1bSJed Brown 
1278c4762a1bSJed Brown TEST*/
1279