xref: /petsc/src/ts/tutorials/ex31.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
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 */
52d71ae5a4SJacob Faibussowitsch PetscErrorCode GetSize(const char *p, PetscInt *sz)
53d71ae5a4SJacob Faibussowitsch {
547510d9b0SBarry Smith   PetscFunctionBeginUser;
5511cc89d2SBarry Smith 
5611cc89d2SBarry Smith   if (!strcmp(p, "hull1972a1") || !strcmp(p, "hull1972a2") || !strcmp(p, "hull1972a3") || !strcmp(p, "hull1972a4") || !strcmp(p, "hull1972a5")) *sz = 1;
5711cc89d2SBarry Smith   else if (!strcmp(p, "hull1972b1")) *sz = 2;
5811cc89d2SBarry Smith   else if (!strcmp(p, "hull1972b2") || !strcmp(p, "hull1972b3") || !strcmp(p, "hull1972b4") || !strcmp(p, "hull1972b5")) *sz = 3;
5911cc89d2SBarry Smith   else if (!strcmp(p, "kulik2013i")) *sz = 4;
6011cc89d2SBarry Smith   else if (!strcmp(p, "hull1972c1") || !strcmp(p, "hull1972c2") || !strcmp(p, "hull1972c3")) *sz = 10;
6111cc89d2SBarry Smith   else if (!strcmp(p, "hull1972c4")) *sz = 52;
6211cc89d2SBarry Smith   else SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_ARG_WRONG, "Incorrect problem name %s", p);
63*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
64c4762a1bSJed Brown }
65c4762a1bSJed Brown 
66c4762a1bSJed Brown /****************************************************************/
67c4762a1bSJed Brown 
68c4762a1bSJed Brown /* Problem specific functions */
69c4762a1bSJed Brown 
70c4762a1bSJed Brown /* Hull, 1972, Problem A1 */
71c4762a1bSJed Brown 
72d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec F, void *s)
73d71ae5a4SJacob Faibussowitsch {
74c4762a1bSJed Brown   PetscScalar       *f;
75c4762a1bSJed Brown   const PetscScalar *y;
76c4762a1bSJed Brown 
777510d9b0SBarry Smith   PetscFunctionBeginUser;
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
80c4762a1bSJed Brown   f[0] = -y[0];
819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
83*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
84c4762a1bSJed Brown }
85c4762a1bSJed Brown 
86d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
87d71ae5a4SJacob Faibussowitsch {
88c4762a1bSJed Brown   const PetscScalar *y;
89c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
90c4762a1bSJed Brown   PetscScalar        value = -1.0;
91c4762a1bSJed Brown 
927510d9b0SBarry Smith   PetscFunctionBeginUser;
939566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
949566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
959566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
969566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
98*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
99c4762a1bSJed Brown }
100c4762a1bSJed Brown 
101d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
102d71ae5a4SJacob Faibussowitsch {
103c4762a1bSJed Brown   const PetscScalar *y;
104c4762a1bSJed Brown   PetscScalar       *f;
105c4762a1bSJed Brown 
1067510d9b0SBarry Smith   PetscFunctionBeginUser;
1079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1089566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
109c4762a1bSJed Brown   f[0] = -y[0];
1109566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
1119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
112c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
1139566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
114*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown 
117d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
118d71ae5a4SJacob Faibussowitsch {
119c4762a1bSJed Brown   const PetscScalar *y;
120c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
121c4762a1bSJed Brown   PetscScalar        value = a - 1.0;
122c4762a1bSJed Brown 
1237510d9b0SBarry Smith   PetscFunctionBeginUser;
1249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1259566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
1269566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1279566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
129*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
130c4762a1bSJed Brown }
131c4762a1bSJed Brown 
132c4762a1bSJed Brown /* Hull, 1972, Problem A2 */
133c4762a1bSJed Brown 
134d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec F, void *s)
135d71ae5a4SJacob Faibussowitsch {
136c4762a1bSJed Brown   const PetscScalar *y;
137c4762a1bSJed Brown   PetscScalar       *f;
138c4762a1bSJed Brown 
1397510d9b0SBarry Smith   PetscFunctionBeginUser;
1409566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1419566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
142c4762a1bSJed Brown   f[0] = -0.5 * y[0] * y[0] * y[0];
1439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
1449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
145*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
146c4762a1bSJed Brown }
147c4762a1bSJed Brown 
148d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
149d71ae5a4SJacob Faibussowitsch {
150c4762a1bSJed Brown   const PetscScalar *y;
151c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
152c4762a1bSJed Brown   PetscScalar        value;
153c4762a1bSJed Brown 
1547510d9b0SBarry Smith   PetscFunctionBeginUser;
1559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
156c4762a1bSJed Brown   value = -0.5 * 3.0 * y[0] * y[0];
1579566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
1589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
161*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
162c4762a1bSJed Brown }
163c4762a1bSJed Brown 
164d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
165d71ae5a4SJacob Faibussowitsch {
166c4762a1bSJed Brown   PetscScalar       *f;
167c4762a1bSJed Brown   const PetscScalar *y;
168c4762a1bSJed Brown 
1697510d9b0SBarry Smith   PetscFunctionBeginUser;
1709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
1719566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
172c4762a1bSJed Brown   f[0] = -0.5 * y[0] * y[0] * y[0];
1739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
1749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
175c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
1769566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
177*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
181d71ae5a4SJacob Faibussowitsch {
182c4762a1bSJed Brown   const PetscScalar *y;
183c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
184c4762a1bSJed Brown   PetscScalar        value;
185c4762a1bSJed Brown 
1867510d9b0SBarry Smith   PetscFunctionBeginUser;
1879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
188c4762a1bSJed Brown   value = a + 0.5 * 3.0 * y[0] * y[0];
1899566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
1909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
193*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown /* Hull, 1972, Problem A3 */
197c4762a1bSJed Brown 
198d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec F, void *s)
199d71ae5a4SJacob Faibussowitsch {
200c4762a1bSJed Brown   const PetscScalar *y;
201c4762a1bSJed Brown   PetscScalar       *f;
202c4762a1bSJed Brown 
2037510d9b0SBarry Smith   PetscFunctionBeginUser;
2049566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2059566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
206c4762a1bSJed Brown   f[0] = y[0] * PetscCosReal(t);
2079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
2089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
209*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
210c4762a1bSJed Brown }
211c4762a1bSJed Brown 
212d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
213d71ae5a4SJacob Faibussowitsch {
214c4762a1bSJed Brown   const PetscScalar *y;
215c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
216c4762a1bSJed Brown   PetscScalar        value = PetscCosReal(t);
217c4762a1bSJed Brown 
2187510d9b0SBarry Smith   PetscFunctionBeginUser;
2199566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2209566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
2219566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2229566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
224*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
225c4762a1bSJed Brown }
226c4762a1bSJed Brown 
227d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
228d71ae5a4SJacob Faibussowitsch {
229c4762a1bSJed Brown   PetscScalar       *f;
230c4762a1bSJed Brown   const PetscScalar *y;
231c4762a1bSJed Brown 
2327510d9b0SBarry Smith   PetscFunctionBeginUser;
2339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2349566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
235c4762a1bSJed Brown   f[0] = y[0] * PetscCosReal(t);
2369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
2379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
238c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
2399566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
240*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
241c4762a1bSJed Brown }
242c4762a1bSJed Brown 
243d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
244d71ae5a4SJacob Faibussowitsch {
245c4762a1bSJed Brown   const PetscScalar *y;
246c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
247c4762a1bSJed Brown   PetscScalar        value = a - PetscCosReal(t);
248c4762a1bSJed Brown 
2497510d9b0SBarry Smith   PetscFunctionBeginUser;
2509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2519566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
2529566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
255*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
256c4762a1bSJed Brown }
257c4762a1bSJed Brown 
258c4762a1bSJed Brown /* Hull, 1972, Problem A4 */
259c4762a1bSJed Brown 
260d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec F, void *s)
261d71ae5a4SJacob Faibussowitsch {
262c4762a1bSJed Brown   PetscScalar       *f;
263c4762a1bSJed Brown   const PetscScalar *y;
264c4762a1bSJed Brown 
2657510d9b0SBarry Smith   PetscFunctionBeginUser;
2669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
268c4762a1bSJed Brown   f[0] = (0.25 * y[0]) * (1.0 - 0.05 * y[0]);
2699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
2709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
271*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
272c4762a1bSJed Brown }
273c4762a1bSJed Brown 
274d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
275d71ae5a4SJacob Faibussowitsch {
276c4762a1bSJed Brown   const PetscScalar *y;
277c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
278c4762a1bSJed Brown   PetscScalar        value;
279c4762a1bSJed Brown 
2807510d9b0SBarry Smith   PetscFunctionBeginUser;
2819566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
282c4762a1bSJed Brown   value = 0.25 * (1.0 - 0.05 * y[0]) - (0.25 * y[0]) * 0.05;
2839566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
2849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
2869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
287*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
288c4762a1bSJed Brown }
289c4762a1bSJed Brown 
290d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
291d71ae5a4SJacob Faibussowitsch {
292c4762a1bSJed Brown   PetscScalar       *f;
293c4762a1bSJed Brown   const PetscScalar *y;
294c4762a1bSJed Brown 
2957510d9b0SBarry Smith   PetscFunctionBeginUser;
2969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
2979566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
298c4762a1bSJed Brown   f[0] = (0.25 * y[0]) * (1.0 - 0.05 * y[0]);
2999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
3009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
301c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
3029566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
303*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
304c4762a1bSJed Brown }
305c4762a1bSJed Brown 
306d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
307d71ae5a4SJacob Faibussowitsch {
308c4762a1bSJed Brown   const PetscScalar *y;
309c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
310c4762a1bSJed Brown   PetscScalar        value;
311c4762a1bSJed Brown 
3127510d9b0SBarry Smith   PetscFunctionBeginUser;
3139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
314c4762a1bSJed Brown   value = a - 0.25 * (1.0 - 0.05 * y[0]) + (0.25 * y[0]) * 0.05;
3159566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
3169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3179566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
319*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
320c4762a1bSJed Brown }
321c4762a1bSJed Brown 
322c4762a1bSJed Brown /* Hull, 1972, Problem A5 */
323c4762a1bSJed Brown 
324d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec F, void *s)
325d71ae5a4SJacob Faibussowitsch {
326c4762a1bSJed Brown   PetscScalar       *f;
327c4762a1bSJed Brown   const PetscScalar *y;
328c4762a1bSJed Brown 
3297510d9b0SBarry Smith   PetscFunctionBeginUser;
3309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
3319566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
332c4762a1bSJed Brown   f[0] = (y[0] - t) / (y[0] + t);
3339566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
3349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
335*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
336c4762a1bSJed Brown }
337c4762a1bSJed Brown 
338d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
339d71ae5a4SJacob Faibussowitsch {
340c4762a1bSJed Brown   const PetscScalar *y;
341c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
342c4762a1bSJed Brown   PetscScalar        value;
343c4762a1bSJed Brown 
3447510d9b0SBarry Smith   PetscFunctionBeginUser;
3459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
346c4762a1bSJed Brown   value = 2 * t / ((t + y[0]) * (t + y[0]));
3479566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
3489566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
351*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
352c4762a1bSJed Brown }
353c4762a1bSJed Brown 
354d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
355d71ae5a4SJacob Faibussowitsch {
356c4762a1bSJed Brown   PetscScalar       *f;
357c4762a1bSJed Brown   const PetscScalar *y;
358c4762a1bSJed Brown 
3597510d9b0SBarry Smith   PetscFunctionBeginUser;
3609566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
3619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
362c4762a1bSJed Brown   f[0] = (y[0] - t) / (y[0] + t);
3639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
3649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
365c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
3669566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
367*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
368c4762a1bSJed Brown }
369c4762a1bSJed Brown 
370d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
371d71ae5a4SJacob Faibussowitsch {
372c4762a1bSJed Brown   const PetscScalar *y;
373c4762a1bSJed Brown   PetscInt           row = 0, col = 0;
374c4762a1bSJed Brown   PetscScalar        value;
375c4762a1bSJed Brown 
3767510d9b0SBarry Smith   PetscFunctionBeginUser;
3779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
378c4762a1bSJed Brown   value = a - 2 * t / ((t + y[0]) * (t + y[0]));
3799566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &row, 1, &col, &value, INSERT_VALUES));
3809566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
3819566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
3829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
383*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
384c4762a1bSJed Brown }
385c4762a1bSJed Brown 
386c4762a1bSJed Brown /* Hull, 1972, Problem B1 */
387c4762a1bSJed Brown 
388d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec F, void *s)
389d71ae5a4SJacob Faibussowitsch {
390c4762a1bSJed Brown   PetscScalar       *f;
391c4762a1bSJed Brown   const PetscScalar *y;
392c4762a1bSJed Brown 
3937510d9b0SBarry Smith   PetscFunctionBeginUser;
3949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
3959566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
396c4762a1bSJed Brown   f[0] = 2.0 * (y[0] - y[0] * y[1]);
397c4762a1bSJed Brown   f[1] = -(y[1] - y[0] * y[1]);
3989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
3999566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
400*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
401c4762a1bSJed Brown }
402c4762a1bSJed Brown 
403d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
404d71ae5a4SJacob Faibussowitsch {
405c4762a1bSJed Brown   PetscScalar       *f;
406c4762a1bSJed Brown   const PetscScalar *y;
407c4762a1bSJed Brown 
4087510d9b0SBarry Smith   PetscFunctionBeginUser;
4099566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4109566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
411c4762a1bSJed Brown   f[0] = 2.0 * (y[0] - y[0] * y[1]);
412c4762a1bSJed Brown   f[1] = -(y[1] - y[0] * y[1]);
4139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
4149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
415c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
4169566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
417*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
418c4762a1bSJed Brown }
419c4762a1bSJed Brown 
420d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
421d71ae5a4SJacob Faibussowitsch {
422c4762a1bSJed Brown   const PetscScalar *y;
423c4762a1bSJed Brown   PetscInt           row[2] = {0, 1};
424c4762a1bSJed Brown   PetscScalar        value[2][2];
425c4762a1bSJed Brown 
4267510d9b0SBarry Smith   PetscFunctionBeginUser;
4279566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4289371c9d4SSatish Balay   value[0][0] = a - 2.0 * (1.0 - y[1]);
4299371c9d4SSatish Balay   value[0][1] = 2.0 * y[0];
4309371c9d4SSatish Balay   value[1][0] = -y[1];
4319371c9d4SSatish Balay   value[1][1] = a + 1.0 - y[0];
4329566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, &row[0], 2, &row[0], &value[0][0], INSERT_VALUES));
4339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4349566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
436*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
437c4762a1bSJed Brown }
438c4762a1bSJed Brown 
439c4762a1bSJed Brown /* Hull, 1972, Problem B2 */
440c4762a1bSJed Brown 
441d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec F, void *s)
442d71ae5a4SJacob Faibussowitsch {
443c4762a1bSJed Brown   PetscScalar       *f;
444c4762a1bSJed Brown   const PetscScalar *y;
445c4762a1bSJed Brown 
4467510d9b0SBarry Smith   PetscFunctionBeginUser;
4479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
449c4762a1bSJed Brown   f[0] = -y[0] + y[1];
450c4762a1bSJed Brown   f[1] = y[0] - 2.0 * y[1] + y[2];
451c4762a1bSJed Brown   f[2] = y[1] - y[2];
4529566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
4539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
454*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
455c4762a1bSJed Brown }
456c4762a1bSJed Brown 
457d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
458d71ae5a4SJacob Faibussowitsch {
459c4762a1bSJed Brown   PetscScalar       *f;
460c4762a1bSJed Brown   const PetscScalar *y;
461c4762a1bSJed Brown 
4627510d9b0SBarry Smith   PetscFunctionBeginUser;
4639566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4649566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
465c4762a1bSJed Brown   f[0] = -y[0] + y[1];
466c4762a1bSJed Brown   f[1] = y[0] - 2.0 * y[1] + y[2];
467c4762a1bSJed Brown   f[2] = y[1] - y[2];
4689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
4699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
470c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
4719566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
472*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
473c4762a1bSJed Brown }
474c4762a1bSJed Brown 
475d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
476d71ae5a4SJacob Faibussowitsch {
477c4762a1bSJed Brown   const PetscScalar *y;
478c4762a1bSJed Brown   PetscInt           row[3] = {0, 1, 2};
479c4762a1bSJed Brown   PetscScalar        value[3][3];
480c4762a1bSJed Brown 
4817510d9b0SBarry Smith   PetscFunctionBeginUser;
4829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
4839371c9d4SSatish Balay   value[0][0] = a + 1.0;
4849371c9d4SSatish Balay   value[0][1] = -1.0;
4859371c9d4SSatish Balay   value[0][2] = 0;
4869371c9d4SSatish Balay   value[1][0] = -1.0;
4879371c9d4SSatish Balay   value[1][1] = a + 2.0;
4889371c9d4SSatish Balay   value[1][2] = -1.0;
4899371c9d4SSatish Balay   value[2][0] = 0;
4909371c9d4SSatish Balay   value[2][1] = -1.0;
4919371c9d4SSatish Balay   value[2][2] = a + 1.0;
4929566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES));
4939566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
4949566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
4959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
496*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
497c4762a1bSJed Brown }
498c4762a1bSJed Brown 
499c4762a1bSJed Brown /* Hull, 1972, Problem B3 */
500c4762a1bSJed Brown 
501d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec F, void *s)
502d71ae5a4SJacob Faibussowitsch {
503c4762a1bSJed Brown   PetscScalar       *f;
504c4762a1bSJed Brown   const PetscScalar *y;
505c4762a1bSJed Brown 
5067510d9b0SBarry Smith   PetscFunctionBeginUser;
5079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
5089566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
509c4762a1bSJed Brown   f[0] = -y[0];
510c4762a1bSJed Brown   f[1] = y[0] - y[1] * y[1];
511c4762a1bSJed Brown   f[2] = y[1] * y[1];
5129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
5139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
514*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
515c4762a1bSJed Brown }
516c4762a1bSJed Brown 
517d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
518d71ae5a4SJacob Faibussowitsch {
519c4762a1bSJed Brown   PetscScalar       *f;
520c4762a1bSJed Brown   const PetscScalar *y;
521c4762a1bSJed Brown 
5227510d9b0SBarry Smith   PetscFunctionBeginUser;
5239566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
5249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
525c4762a1bSJed Brown   f[0] = -y[0];
526c4762a1bSJed Brown   f[1] = y[0] - y[1] * y[1];
527c4762a1bSJed Brown   f[2] = y[1] * y[1];
5289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
5299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
530c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
5319566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
532*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
533c4762a1bSJed Brown }
534c4762a1bSJed Brown 
535d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
536d71ae5a4SJacob Faibussowitsch {
537c4762a1bSJed Brown   const PetscScalar *y;
538c4762a1bSJed Brown   PetscInt           row[3] = {0, 1, 2};
539c4762a1bSJed Brown   PetscScalar        value[3][3];
540c4762a1bSJed Brown 
5417510d9b0SBarry Smith   PetscFunctionBeginUser;
5429566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
5439371c9d4SSatish Balay   value[0][0] = a + 1.0;
5449371c9d4SSatish Balay   value[0][1] = 0;
5459371c9d4SSatish Balay   value[0][2] = 0;
5469371c9d4SSatish Balay   value[1][0] = -1.0;
5479371c9d4SSatish Balay   value[1][1] = a + 2.0 * y[1];
5489371c9d4SSatish Balay   value[1][2] = 0;
5499371c9d4SSatish Balay   value[2][0] = 0;
5509371c9d4SSatish Balay   value[2][1] = -2.0 * y[1];
5519371c9d4SSatish Balay   value[2][2] = a;
5529566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES));
5539566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
5549566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
5559566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
556*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
557c4762a1bSJed Brown }
558c4762a1bSJed Brown 
559c4762a1bSJed Brown /* Hull, 1972, Problem B4 */
560c4762a1bSJed Brown 
561d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec F, void *s)
562d71ae5a4SJacob Faibussowitsch {
563c4762a1bSJed Brown   PetscScalar       *f;
564c4762a1bSJed Brown   const PetscScalar *y;
565c4762a1bSJed Brown 
5667510d9b0SBarry Smith   PetscFunctionBeginUser;
5679566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
5689566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
569c4762a1bSJed Brown   f[0] = -y[1] - y[0] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
570c4762a1bSJed Brown   f[1] = y[0] - y[1] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
571c4762a1bSJed Brown   f[2] = y[0] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
5729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
5739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
574*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
575c4762a1bSJed Brown }
576c4762a1bSJed Brown 
577d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
578d71ae5a4SJacob Faibussowitsch {
579c4762a1bSJed Brown   PetscScalar       *f;
580c4762a1bSJed Brown   const PetscScalar *y;
581c4762a1bSJed Brown 
5827510d9b0SBarry Smith   PetscFunctionBeginUser;
5839566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
5849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
585c4762a1bSJed Brown   f[0] = -y[1] - y[0] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
586c4762a1bSJed Brown   f[1] = y[0] - y[1] * y[2] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
587c4762a1bSJed Brown   f[2] = y[0] / PetscSqrtScalar(y[0] * y[0] + y[1] * y[1]);
5889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
5899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
590c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
5919566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
592*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
593c4762a1bSJed Brown }
594c4762a1bSJed Brown 
595d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
596d71ae5a4SJacob Faibussowitsch {
597c4762a1bSJed Brown   const PetscScalar *y;
598c4762a1bSJed Brown   PetscInt           row[3] = {0, 1, 2};
599c4762a1bSJed Brown   PetscScalar        value[3][3], fac, fac2;
600c4762a1bSJed Brown 
6017510d9b0SBarry Smith   PetscFunctionBeginUser;
6029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
603c4762a1bSJed Brown   fac         = PetscPowScalar(y[0] * y[0] + y[1] * y[1], -1.5);
604c4762a1bSJed Brown   fac2        = PetscPowScalar(y[0] * y[0] + y[1] * y[1], -0.5);
605c4762a1bSJed Brown   value[0][0] = a + (y[1] * y[1] * y[2]) * fac;
606c4762a1bSJed Brown   value[0][1] = 1.0 - (y[0] * y[1] * y[2]) * fac;
607c4762a1bSJed Brown   value[0][2] = y[0] * fac2;
608c4762a1bSJed Brown   value[1][0] = -1.0 - y[0] * y[1] * y[2] * fac;
609c4762a1bSJed Brown   value[1][1] = a + y[0] * y[0] * y[2] * fac;
610c4762a1bSJed Brown   value[1][2] = y[1] * fac2;
611c4762a1bSJed Brown   value[2][0] = -y[1] * y[1] * fac;
612c4762a1bSJed Brown   value[2][1] = y[0] * y[1] * fac;
613c4762a1bSJed Brown   value[2][2] = a;
6149566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES));
6159566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
6169566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
6179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
618*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
619c4762a1bSJed Brown }
620c4762a1bSJed Brown 
621c4762a1bSJed Brown /* Hull, 1972, Problem B5 */
622c4762a1bSJed Brown 
623d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec F, void *s)
624d71ae5a4SJacob Faibussowitsch {
625c4762a1bSJed Brown   PetscScalar       *f;
626c4762a1bSJed Brown   const PetscScalar *y;
627c4762a1bSJed Brown 
6287510d9b0SBarry Smith   PetscFunctionBeginUser;
6299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
6309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
631c4762a1bSJed Brown   f[0] = y[1] * y[2];
632c4762a1bSJed Brown   f[1] = -y[0] * y[2];
633c4762a1bSJed Brown   f[2] = -0.51 * y[0] * y[1];
6349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
6359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
636*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
637c4762a1bSJed Brown }
638c4762a1bSJed Brown 
639d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
640d71ae5a4SJacob Faibussowitsch {
641c4762a1bSJed Brown   PetscScalar       *f;
642c4762a1bSJed Brown   const PetscScalar *y;
643c4762a1bSJed Brown 
6447510d9b0SBarry Smith   PetscFunctionBeginUser;
6459566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
6469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
647c4762a1bSJed Brown   f[0] = y[1] * y[2];
648c4762a1bSJed Brown   f[1] = -y[0] * y[2];
649c4762a1bSJed Brown   f[2] = -0.51 * y[0] * y[1];
6509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
6519566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
652c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
6539566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
654*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
655c4762a1bSJed Brown }
656c4762a1bSJed Brown 
657d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
658d71ae5a4SJacob Faibussowitsch {
659c4762a1bSJed Brown   const PetscScalar *y;
660c4762a1bSJed Brown   PetscInt           row[3] = {0, 1, 2};
661c4762a1bSJed Brown   PetscScalar        value[3][3];
662c4762a1bSJed Brown 
6637510d9b0SBarry Smith   PetscFunctionBeginUser;
6649566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
6659371c9d4SSatish Balay   value[0][0] = a;
6669371c9d4SSatish Balay   value[0][1] = -y[2];
6679371c9d4SSatish Balay   value[0][2] = -y[1];
6689371c9d4SSatish Balay   value[1][0] = y[2];
6699371c9d4SSatish Balay   value[1][1] = a;
6709371c9d4SSatish Balay   value[1][2] = y[0];
6719371c9d4SSatish Balay   value[2][0] = 0.51 * y[1];
6729371c9d4SSatish Balay   value[2][1] = 0.51 * y[0];
6739371c9d4SSatish Balay   value[2][2] = a;
6749566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 3, &row[0], 3, &row[0], &value[0][0], INSERT_VALUES));
6759566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
6769566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
6779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
678*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
679c4762a1bSJed Brown }
680c4762a1bSJed Brown 
681c4762a1bSJed Brown /* Kulikov, 2013, Problem I */
682c4762a1bSJed Brown 
683d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec F, void *s)
684d71ae5a4SJacob Faibussowitsch {
685c4762a1bSJed Brown   PetscScalar       *f;
686c4762a1bSJed Brown   const PetscScalar *y;
687c4762a1bSJed Brown 
6887510d9b0SBarry Smith   PetscFunctionBeginUser;
6899566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
6909566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
691c4762a1bSJed Brown   f[0] = 2. * t * PetscPowScalar(y[1], 1. / 5.) * y[3];
692c4762a1bSJed Brown   f[1] = 10. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.));
693c4762a1bSJed Brown   f[2] = 2. * t * y[3];
694c4762a1bSJed Brown   f[3] = -2. * t * PetscLogScalar(y[0]);
6959566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
6969566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
697*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
698c4762a1bSJed Brown }
699c4762a1bSJed Brown 
700d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
701d71ae5a4SJacob Faibussowitsch {
702c4762a1bSJed Brown   const PetscScalar *y;
703c4762a1bSJed Brown   PetscInt           row[4] = {0, 1, 2, 3};
704c4762a1bSJed Brown   PetscScalar        value[4][4];
705c4762a1bSJed Brown   PetscScalar        m1, m2;
7067510d9b0SBarry Smith   PetscFunctionBeginUser;
7079566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
708c4762a1bSJed Brown   m1          = (2. * t * y[3]) / (5. * PetscPowScalar(y[1], 4. / 5.));
709c4762a1bSJed Brown   m2          = 2. * t * PetscPowScalar(y[1], 1. / 5.);
7109371c9d4SSatish Balay   value[0][0] = 0.;
7119371c9d4SSatish Balay   value[0][1] = m1;
7129371c9d4SSatish Balay   value[0][2] = 0.;
7139371c9d4SSatish Balay   value[0][3] = m2;
714c4762a1bSJed Brown   m1          = 50. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.));
715c4762a1bSJed Brown   m2          = 10. * t * PetscExpScalar(5.0 * (y[2] - 1.));
7169371c9d4SSatish Balay   value[1][0] = 0.;
7179371c9d4SSatish Balay   value[1][1] = 0.;
7189371c9d4SSatish Balay   value[1][2] = m1;
7199371c9d4SSatish Balay   value[1][3] = m2;
7209371c9d4SSatish Balay   value[2][0] = 0.;
7219371c9d4SSatish Balay   value[2][1] = 0.;
7229371c9d4SSatish Balay   value[2][2] = 0.;
7239371c9d4SSatish Balay   value[2][3] = 2 * t;
7249371c9d4SSatish Balay   value[3][0] = -2. * t / y[0];
7259371c9d4SSatish Balay   value[3][1] = 0.;
7269371c9d4SSatish Balay   value[3][2] = 0.;
7279371c9d4SSatish Balay   value[3][3] = 0.;
7289566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 4, &row[0], 4, &row[0], &value[0][0], INSERT_VALUES));
7299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
7319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
732*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
733c4762a1bSJed Brown }
734c4762a1bSJed Brown 
735d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
736d71ae5a4SJacob Faibussowitsch {
737c4762a1bSJed Brown   PetscScalar       *f;
738c4762a1bSJed Brown   const PetscScalar *y;
739c4762a1bSJed Brown 
7407510d9b0SBarry Smith   PetscFunctionBeginUser;
7419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
7429566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
743c4762a1bSJed Brown   f[0] = 2. * t * PetscPowScalar(y[1], 1. / 5.) * y[3];
744c4762a1bSJed Brown   f[1] = 10. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.));
745c4762a1bSJed Brown   f[2] = 2. * t * y[3];
746c4762a1bSJed Brown   f[3] = -2. * t * PetscLogScalar(y[0]);
7479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
7489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
749c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
7509566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
751*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
752c4762a1bSJed Brown }
753c4762a1bSJed Brown 
754d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
755d71ae5a4SJacob Faibussowitsch {
756c4762a1bSJed Brown   const PetscScalar *y;
757c4762a1bSJed Brown   PetscInt           row[4] = {0, 1, 2, 3};
758c4762a1bSJed Brown   PetscScalar        value[4][4];
759c4762a1bSJed Brown   PetscScalar        m1, m2;
760c4762a1bSJed Brown 
7617510d9b0SBarry Smith   PetscFunctionBeginUser;
7629566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
763c4762a1bSJed Brown   m1          = (2. * t * y[3]) / (5. * PetscPowScalar(y[1], 4. / 5.));
764c4762a1bSJed Brown   m2          = 2. * t * PetscPowScalar(y[1], 1. / 5.);
7659371c9d4SSatish Balay   value[0][0] = a;
7669371c9d4SSatish Balay   value[0][1] = m1;
7679371c9d4SSatish Balay   value[0][2] = 0.;
7689371c9d4SSatish Balay   value[0][3] = m2;
769c4762a1bSJed Brown   m1          = 50. * t * y[3] * PetscExpScalar(5.0 * (y[2] - 1.));
770c4762a1bSJed Brown   m2          = 10. * t * PetscExpScalar(5.0 * (y[2] - 1.));
7719371c9d4SSatish Balay   value[1][0] = 0.;
7729371c9d4SSatish Balay   value[1][1] = a;
7739371c9d4SSatish Balay   value[1][2] = m1;
7749371c9d4SSatish Balay   value[1][3] = m2;
7759371c9d4SSatish Balay   value[2][0] = 0.;
7769371c9d4SSatish Balay   value[2][1] = 0.;
7779371c9d4SSatish Balay   value[2][2] = a;
7789371c9d4SSatish Balay   value[2][3] = 2 * t;
7799371c9d4SSatish Balay   value[3][0] = -2. * t / y[0];
7809371c9d4SSatish Balay   value[3][1] = 0.;
7819371c9d4SSatish Balay   value[3][2] = 0.;
7829371c9d4SSatish Balay   value[3][3] = a;
7839566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 4, &row[0], 4, &row[0], &value[0][0], INSERT_VALUES));
7849566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
7859566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
7869566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
787*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
788c4762a1bSJed Brown }
789c4762a1bSJed Brown 
790c4762a1bSJed Brown /* Hull, 1972, Problem C1 */
791c4762a1bSJed Brown 
792d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec F, void *s)
793d71ae5a4SJacob Faibussowitsch {
794c4762a1bSJed Brown   PetscScalar       *f;
795c4762a1bSJed Brown   const PetscScalar *y;
796c4762a1bSJed Brown   PetscInt           N, i;
797c4762a1bSJed Brown 
7987510d9b0SBarry Smith   PetscFunctionBeginUser;
7999566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8009566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
8019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
802c4762a1bSJed Brown   f[0] = -y[0];
803ad540459SPierre Jolivet   for (i = 1; i < N - 1; i++) f[i] = y[i - 1] - y[i];
804c4762a1bSJed Brown   f[N - 1] = y[N - 2];
8059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
8069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
807*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
808c4762a1bSJed Brown }
809c4762a1bSJed Brown 
810d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
811d71ae5a4SJacob Faibussowitsch {
812c4762a1bSJed Brown   PetscScalar       *f;
813c4762a1bSJed Brown   const PetscScalar *y;
814c4762a1bSJed Brown   PetscInt           N, i;
815c4762a1bSJed Brown 
8167510d9b0SBarry Smith   PetscFunctionBeginUser;
8179566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
8199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
820c4762a1bSJed Brown   f[0] = -y[0];
821ad540459SPierre Jolivet   for (i = 1; i < N - 1; i++) f[i] = y[i - 1] - y[i];
822c4762a1bSJed Brown   f[N - 1] = y[N - 2];
8239566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
8249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
825c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
8269566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
827*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
828c4762a1bSJed Brown }
829c4762a1bSJed Brown 
830d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
831d71ae5a4SJacob Faibussowitsch {
832c4762a1bSJed Brown   const PetscScalar *y;
833c4762a1bSJed Brown   PetscInt           N, i, col[2];
834c4762a1bSJed Brown   PetscScalar        value[2];
835c4762a1bSJed Brown 
8367510d9b0SBarry Smith   PetscFunctionBeginUser;
8379566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
839c4762a1bSJed Brown   i        = 0;
8409371c9d4SSatish Balay   value[0] = a + 1;
8419371c9d4SSatish Balay   col[0]   = 0;
8429371c9d4SSatish Balay   value[1] = 0;
8439371c9d4SSatish Balay   col[1]   = 1;
8449566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
845c4762a1bSJed Brown   for (i = 0; i < N; i++) {
8469371c9d4SSatish Balay     value[0] = -1;
8479371c9d4SSatish Balay     col[0]   = i - 1;
8489371c9d4SSatish Balay     value[1] = a + 1;
8499371c9d4SSatish Balay     col[1]   = i;
8509566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
851c4762a1bSJed Brown   }
852c4762a1bSJed Brown   i        = N - 1;
8539371c9d4SSatish Balay   value[0] = -1;
8549371c9d4SSatish Balay   col[0]   = N - 2;
8559371c9d4SSatish Balay   value[1] = a;
8569371c9d4SSatish Balay   col[1]   = N - 1;
8579566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
8589566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
8599566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
8609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
861*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
862c4762a1bSJed Brown }
863c4762a1bSJed Brown 
864c4762a1bSJed Brown /* Hull, 1972, Problem C2 */
865c4762a1bSJed Brown 
866d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec F, void *s)
867d71ae5a4SJacob Faibussowitsch {
868c4762a1bSJed Brown   const PetscScalar *y;
869c4762a1bSJed Brown   PetscScalar       *f;
870c4762a1bSJed Brown   PetscInt           N, i;
871c4762a1bSJed Brown 
8727510d9b0SBarry Smith   PetscFunctionBeginUser;
8739566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8749566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
8759566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
876c4762a1bSJed Brown   f[0] = -y[0];
877ad540459SPierre Jolivet   for (i = 1; i < N - 1; i++) f[i] = (PetscReal)i * y[i - 1] - (PetscReal)(i + 1) * y[i];
878c4762a1bSJed Brown   f[N - 1] = (PetscReal)(N - 1) * y[N - 2];
8799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
8809566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
881*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
882c4762a1bSJed Brown }
883c4762a1bSJed Brown 
884d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
885d71ae5a4SJacob Faibussowitsch {
886c4762a1bSJed Brown   PetscScalar       *f;
887c4762a1bSJed Brown   const PetscScalar *y;
888c4762a1bSJed Brown   PetscInt           N, i;
889c4762a1bSJed Brown 
8907510d9b0SBarry Smith   PetscFunctionBeginUser;
8919566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
8929566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
8939566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
894c4762a1bSJed Brown   f[0] = -y[0];
895ad540459SPierre Jolivet   for (i = 1; i < N - 1; i++) f[i] = (PetscReal)i * y[i - 1] - (PetscReal)(i + 1) * y[i];
896c4762a1bSJed Brown   f[N - 1] = (PetscReal)(N - 1) * y[N - 2];
8979566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
8989566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
899c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
9009566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
901*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
902c4762a1bSJed Brown }
903c4762a1bSJed Brown 
904d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
905d71ae5a4SJacob Faibussowitsch {
906c4762a1bSJed Brown   const PetscScalar *y;
907c4762a1bSJed Brown   PetscInt           N, i, col[2];
908c4762a1bSJed Brown   PetscScalar        value[2];
909c4762a1bSJed Brown 
9107510d9b0SBarry Smith   PetscFunctionBeginUser;
9119566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
9129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
913c4762a1bSJed Brown   i        = 0;
9149371c9d4SSatish Balay   value[0] = a + 1;
9159371c9d4SSatish Balay   col[0]   = 0;
9169371c9d4SSatish Balay   value[1] = 0;
9179371c9d4SSatish Balay   col[1]   = 1;
9189566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
919c4762a1bSJed Brown   for (i = 0; i < N; i++) {
9209371c9d4SSatish Balay     value[0] = -(PetscReal)i;
9219371c9d4SSatish Balay     col[0]   = i - 1;
9229371c9d4SSatish Balay     value[1] = a + (PetscReal)(i + 1);
9239371c9d4SSatish Balay     col[1]   = i;
9249566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
925c4762a1bSJed Brown   }
926c4762a1bSJed Brown   i        = N - 1;
9279371c9d4SSatish Balay   value[0] = -(PetscReal)(N - 1);
9289371c9d4SSatish Balay   col[0]   = N - 2;
9299371c9d4SSatish Balay   value[1] = a;
9309371c9d4SSatish Balay   col[1]   = N - 1;
9319566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 1, &i, 2, col, value, INSERT_VALUES));
9329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
9339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
9349566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
935*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
936c4762a1bSJed Brown }
937c4762a1bSJed Brown 
938c4762a1bSJed Brown /* Hull, 1972, Problem C3 and C4 */
939c4762a1bSJed Brown 
940d71ae5a4SJacob Faibussowitsch PetscErrorCode RHSFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec F, void *s)
941d71ae5a4SJacob Faibussowitsch {
942c4762a1bSJed Brown   PetscScalar       *f;
943c4762a1bSJed Brown   const PetscScalar *y;
944c4762a1bSJed Brown   PetscInt           N, i;
945c4762a1bSJed Brown 
9467510d9b0SBarry Smith   PetscFunctionBeginUser;
9479566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
9489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
9499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
950c4762a1bSJed Brown   f[0] = -2.0 * y[0] + y[1];
951ad540459SPierre Jolivet   for (i = 1; i < N - 1; i++) f[i] = y[i - 1] - 2.0 * y[i] + y[i + 1];
952c4762a1bSJed Brown   f[N - 1] = y[N - 2] - 2.0 * y[N - 1];
9539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
9549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
955*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
956c4762a1bSJed Brown }
957c4762a1bSJed Brown 
958d71ae5a4SJacob Faibussowitsch PetscErrorCode IFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
959d71ae5a4SJacob Faibussowitsch {
960c4762a1bSJed Brown   PetscScalar       *f;
961c4762a1bSJed Brown   const PetscScalar *y;
962c4762a1bSJed Brown   PetscInt           N, i;
963c4762a1bSJed Brown 
9647510d9b0SBarry Smith   PetscFunctionBeginUser;
9659566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
9669566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
9679566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
968c4762a1bSJed Brown   f[0] = -2.0 * y[0] + y[1];
969ad540459SPierre Jolivet   for (i = 1; i < N - 1; i++) f[i] = y[i - 1] - 2.0 * y[i] + y[i + 1];
970c4762a1bSJed Brown   f[N - 1] = y[N - 2] - 2.0 * y[N - 1];
9719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
9729566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
973c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
9749566063dSJacob Faibussowitsch   PetscCall(VecAYPX(F, -1.0, Ydot));
975*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
976c4762a1bSJed Brown }
977c4762a1bSJed Brown 
978d71ae5a4SJacob Faibussowitsch PetscErrorCode IJacobian_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
979d71ae5a4SJacob Faibussowitsch {
980c4762a1bSJed Brown   const PetscScalar *y;
981c4762a1bSJed Brown   PetscScalar        value[3];
982c4762a1bSJed Brown   PetscInt           N, i, col[3];
983c4762a1bSJed Brown 
9847510d9b0SBarry Smith   PetscFunctionBeginUser;
9859566063dSJacob Faibussowitsch   PetscCall(VecGetSize(Y, &N));
9869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &y));
987c4762a1bSJed Brown   for (i = 0; i < N; i++) {
988c4762a1bSJed Brown     if (i == 0) {
9899371c9d4SSatish Balay       value[0] = a + 2;
9909371c9d4SSatish Balay       col[0]   = i;
9919371c9d4SSatish Balay       value[1] = -1;
9929371c9d4SSatish Balay       col[1]   = i + 1;
9939371c9d4SSatish Balay       value[2] = 0;
9949371c9d4SSatish Balay       col[2]   = i + 2;
995c4762a1bSJed Brown     } else if (i == N - 1) {
9969371c9d4SSatish Balay       value[0] = 0;
9979371c9d4SSatish Balay       col[0]   = i - 2;
9989371c9d4SSatish Balay       value[1] = -1;
9999371c9d4SSatish Balay       col[1]   = i - 1;
10009371c9d4SSatish Balay       value[2] = a + 2;
10019371c9d4SSatish Balay       col[2]   = i;
1002c4762a1bSJed Brown     } else {
10039371c9d4SSatish Balay       value[0] = -1;
10049371c9d4SSatish Balay       col[0]   = i - 1;
10059371c9d4SSatish Balay       value[1] = a + 2;
10069371c9d4SSatish Balay       col[1]   = i;
10079371c9d4SSatish Balay       value[2] = -1;
10089371c9d4SSatish Balay       col[2]   = i + 1;
1009c4762a1bSJed Brown     }
10109566063dSJacob Faibussowitsch     PetscCall(MatSetValues(A, 1, &i, 3, col, value, INSERT_VALUES));
1011c4762a1bSJed Brown   }
10129566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
10139566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
10149566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &y));
1015*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1016c4762a1bSJed Brown }
1017c4762a1bSJed Brown 
1018c4762a1bSJed Brown /***************************************************************************/
1019c4762a1bSJed Brown 
1020c4762a1bSJed Brown /* Sets the initial solution for the IVP and sets up the function pointers*/
1021d71ae5a4SJacob Faibussowitsch PetscErrorCode Initialize(Vec Y, void *s)
1022d71ae5a4SJacob Faibussowitsch {
1023c4762a1bSJed Brown   char        *p = (char *)s;
1024c4762a1bSJed Brown   PetscScalar *y;
1025c4762a1bSJed Brown   PetscReal    t0;
102611cc89d2SBarry Smith   PetscInt     N;
1027c4762a1bSJed Brown   PetscBool    flg;
1028c4762a1bSJed Brown 
10297510d9b0SBarry Smith   PetscFunctionBeginUser;
103011cc89d2SBarry Smith   PetscCall(GetSize((const char *)s, &N));
1031*3ba16761SJacob Faibussowitsch   PetscCall(VecZeroEntries(Y));
10329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(Y, &y));
1033c4762a1bSJed Brown   if (!strcmp(p, "hull1972a1")) {
1034c4762a1bSJed Brown     y[0]        = 1.0;
1035c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A1;
1036c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A1;
1037c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A1;
1038c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A1;
1039c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a2")) {
1040c4762a1bSJed Brown     y[0]        = 1.0;
1041c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A2;
1042c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A2;
1043c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A2;
1044c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A2;
1045c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a3")) {
1046c4762a1bSJed Brown     y[0]        = 1.0;
1047c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A3;
1048c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A3;
1049c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A3;
1050c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A3;
1051c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a4")) {
1052c4762a1bSJed Brown     y[0]        = 1.0;
1053c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A4;
1054c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A4;
1055c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A4;
1056c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A4;
1057c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a5")) {
1058c4762a1bSJed Brown     y[0]        = 4.0;
1059c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A5;
1060c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A5;
1061c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A5;
1062c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A5;
1063c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b1")) {
1064c4762a1bSJed Brown     y[0]        = 1.0;
1065c4762a1bSJed Brown     y[1]        = 3.0;
1066c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B1;
1067c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B1;
1068c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B1;
1069c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b2")) {
1070c4762a1bSJed Brown     y[0]        = 2.0;
1071c4762a1bSJed Brown     y[1]        = 0.0;
1072c4762a1bSJed Brown     y[2]        = 1.0;
1073c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B2;
1074c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B2;
1075c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B2;
1076c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b3")) {
1077c4762a1bSJed Brown     y[0]        = 1.0;
1078c4762a1bSJed Brown     y[1]        = 0.0;
1079c4762a1bSJed Brown     y[2]        = 0.0;
1080c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B3;
1081c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B3;
1082c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B3;
1083c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b4")) {
1084c4762a1bSJed Brown     y[0]        = 3.0;
1085c4762a1bSJed Brown     y[1]        = 0.0;
1086c4762a1bSJed Brown     y[2]        = 0.0;
1087c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B4;
1088c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B4;
1089c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B4;
1090c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972b5")) {
1091c4762a1bSJed Brown     y[0]        = 0.0;
1092c4762a1bSJed Brown     y[1]        = 1.0;
1093c4762a1bSJed Brown     y[2]        = 1.0;
1094c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B5;
1095c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B5;
1096c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B5;
1097c4762a1bSJed Brown   } else if (!strcmp(p, "kulik2013i")) {
1098c4762a1bSJed Brown     t0          = 0.;
1099c4762a1bSJed Brown     y[0]        = PetscExpReal(PetscSinReal(t0 * t0));
1100c4762a1bSJed Brown     y[1]        = PetscExpReal(5. * PetscSinReal(t0 * t0));
1101c4762a1bSJed Brown     y[2]        = PetscSinReal(t0 * t0) + 1.0;
1102c4762a1bSJed Brown     y[3]        = PetscCosReal(t0 * t0);
1103c4762a1bSJed Brown     RHSFunction = RHSFunction_Kulikov2013I;
1104c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Kulikov2013I;
1105c4762a1bSJed Brown     IFunction   = IFunction_Kulikov2013I;
1106c4762a1bSJed Brown     IJacobian   = IJacobian_Kulikov2013I;
1107c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972c1")) {
1108c4762a1bSJed Brown     y[0]        = 1.0;
1109c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C1;
1110c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C1;
1111c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C1;
1112c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972c2")) {
1113c4762a1bSJed Brown     y[0]        = 1.0;
1114c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C2;
1115c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C2;
1116c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C2;
111711cc89d2SBarry Smith   } else if (!strcmp(p, "hull1972c3") || !strcmp(p, "hull1972c4")) {
1118c4762a1bSJed Brown     y[0]        = 1.0;
1119c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C34;
1120c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C34;
1121c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C34;
1122c4762a1bSJed Brown   }
11239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetScalarArray(NULL, NULL, "-yinit", y, &N, &flg));
11249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(Y, &y));
1125*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1126c4762a1bSJed Brown }
1127c4762a1bSJed Brown 
1128c4762a1bSJed Brown /* Calculates the exact solution to problems that have one */
1129d71ae5a4SJacob Faibussowitsch PetscErrorCode ExactSolution(Vec Y, void *s, PetscReal t, PetscBool *flag)
1130d71ae5a4SJacob Faibussowitsch {
1131c4762a1bSJed Brown   char        *p = (char *)s;
1132c4762a1bSJed Brown   PetscScalar *y;
1133c4762a1bSJed Brown 
11347510d9b0SBarry Smith   PetscFunctionBeginUser;
1135c4762a1bSJed Brown   if (!strcmp(p, "hull1972a1")) {
11369566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1137c4762a1bSJed Brown     y[0]  = PetscExpReal(-t);
1138c4762a1bSJed Brown     *flag = PETSC_TRUE;
11399566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1140c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a2")) {
11419566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1142c4762a1bSJed Brown     y[0]  = 1.0 / PetscSqrtReal(t + 1);
1143c4762a1bSJed Brown     *flag = PETSC_TRUE;
11449566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1145c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a3")) {
11469566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1147c4762a1bSJed Brown     y[0]  = PetscExpReal(PetscSinReal(t));
1148c4762a1bSJed Brown     *flag = PETSC_TRUE;
11499566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1150c4762a1bSJed Brown   } else if (!strcmp(p, "hull1972a4")) {
11519566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1152c4762a1bSJed Brown     y[0]  = 20.0 / (1 + 19.0 * PetscExpReal(-t / 4.0));
1153c4762a1bSJed Brown     *flag = PETSC_TRUE;
11549566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1155c4762a1bSJed Brown   } else if (!strcmp(p, "kulik2013i")) {
11569566063dSJacob Faibussowitsch     PetscCall(VecGetArray(Y, &y));
1157c4762a1bSJed Brown     y[0]  = PetscExpReal(PetscSinReal(t * t));
1158c4762a1bSJed Brown     y[1]  = PetscExpReal(5. * PetscSinReal(t * t));
1159c4762a1bSJed Brown     y[2]  = PetscSinReal(t * t) + 1.0;
1160c4762a1bSJed Brown     y[3]  = PetscCosReal(t * t);
1161c4762a1bSJed Brown     *flag = PETSC_TRUE;
11629566063dSJacob Faibussowitsch     PetscCall(VecRestoreArray(Y, &y));
1163c4762a1bSJed Brown   } else {
11649566063dSJacob Faibussowitsch     PetscCall(VecSet(Y, 0));
1165c4762a1bSJed Brown     *flag = PETSC_FALSE;
1166c4762a1bSJed Brown   }
1167*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1168c4762a1bSJed Brown }
1169c4762a1bSJed Brown 
1170c4762a1bSJed Brown /* Solves the specified ODE and computes the error if exact solution is available */
1171d71ae5a4SJacob Faibussowitsch PetscErrorCode SolveODE(char *ptype, PetscReal dt, PetscReal tfinal, PetscInt maxiter, PetscReal *error, PetscBool *exact_flag)
1172d71ae5a4SJacob Faibussowitsch {
1173c4762a1bSJed Brown   TS        ts;          /* time-integrator                        */
1174c4762a1bSJed Brown   Vec       Y;           /* Solution vector                        */
1175c4762a1bSJed Brown   Vec       Yex;         /* Exact solution                         */
1176c4762a1bSJed Brown   PetscInt  N;           /* Size of the system of equations        */
1177c4762a1bSJed Brown   TSType    time_scheme; /* Type of time-integration scheme        */
1178c4762a1bSJed Brown   Mat       Jac = NULL;  /* Jacobian matrix                        */
1179c4762a1bSJed Brown   Vec       Yerr;        /* Auxiliary solution vector              */
1180c4762a1bSJed Brown   PetscReal err_norm;    /* Estimated error norm                   */
1181c4762a1bSJed Brown   PetscReal final_time;  /* Actual final time from the integrator  */
1182c4762a1bSJed Brown 
11837510d9b0SBarry Smith   PetscFunctionBeginUser;
118411cc89d2SBarry Smith   PetscCall(GetSize((const char *)&ptype[0], &N));
11853c633725SBarry Smith   PetscCheck(N >= 0, PETSC_COMM_WORLD, PETSC_ERR_ARG_SIZ, "Illegal problem specification.");
11869566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_WORLD, &Y));
11879566063dSJacob Faibussowitsch   PetscCall(VecSetSizes(Y, N, PETSC_DECIDE));
11889566063dSJacob Faibussowitsch   PetscCall(VecSetUp(Y));
11899566063dSJacob Faibussowitsch   PetscCall(VecSet(Y, 0));
1190c4762a1bSJed Brown 
1191c4762a1bSJed Brown   /* Initialize the problem */
11929566063dSJacob Faibussowitsch   PetscCall(Initialize(Y, &ptype[0]));
1193c4762a1bSJed Brown 
1194c4762a1bSJed Brown   /* Create and initialize the time-integrator                            */
11959566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1196c4762a1bSJed Brown   /* Default time integration options                                     */
11979566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSRK));
11989566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, maxiter));
11999566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, tfinal));
12009566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, dt));
12019566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1202c4762a1bSJed Brown   /* Read command line options for time integration                       */
12039566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1204c4762a1bSJed Brown   /* Set solution vector                                                  */
12059566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, Y));
1206c4762a1bSJed Brown   /* Specify left/right-hand side functions                               */
12079566063dSJacob Faibussowitsch   PetscCall(TSGetType(ts, &time_scheme));
1208c4762a1bSJed Brown 
1209c4762a1bSJed Brown   if ((!strcmp(time_scheme, TSEULER)) || (!strcmp(time_scheme, TSRK)) || (!strcmp(time_scheme, TSSSP) || (!strcmp(time_scheme, TSGLEE)))) {
1210c4762a1bSJed Brown     /* Explicit time-integration -> specify right-hand side function ydot = f(y) */
12119566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &ptype[0]));
12129566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &Jac));
12139566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Jac, PETSC_DECIDE, PETSC_DECIDE, N, N));
12149566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(Jac));
12159566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Jac));
12169566063dSJacob Faibussowitsch     PetscCall(TSSetRHSJacobian(ts, Jac, Jac, RHSJacobian, &ptype[0]));
1217c4762a1bSJed Brown   } else if ((!strcmp(time_scheme, TSTHETA)) || (!strcmp(time_scheme, TSBEULER)) || (!strcmp(time_scheme, TSCN)) || (!strcmp(time_scheme, TSALPHA)) || (!strcmp(time_scheme, TSARKIMEX))) {
1218c4762a1bSJed Brown     /* Implicit time-integration -> specify left-hand side function ydot-f(y) = 0 */
1219c4762a1bSJed Brown     /* and its Jacobian function                                                 */
12209566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction, &ptype[0]));
12219566063dSJacob Faibussowitsch     PetscCall(MatCreate(PETSC_COMM_WORLD, &Jac));
12229566063dSJacob Faibussowitsch     PetscCall(MatSetSizes(Jac, PETSC_DECIDE, PETSC_DECIDE, N, N));
12239566063dSJacob Faibussowitsch     PetscCall(MatSetFromOptions(Jac));
12249566063dSJacob Faibussowitsch     PetscCall(MatSetUp(Jac));
12259566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, Jac, Jac, IJacobian, &ptype[0]));
1226c4762a1bSJed Brown   }
1227c4762a1bSJed Brown 
1228c4762a1bSJed Brown   /* Solve */
12299566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, Y));
12309566063dSJacob Faibussowitsch   PetscCall(TSGetTime(ts, &final_time));
1231c4762a1bSJed Brown 
1232c4762a1bSJed Brown   /* Get the estimated error, if available */
12339566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(Y, &Yerr));
12349566063dSJacob Faibussowitsch   PetscCall(VecZeroEntries(Yerr));
12359566063dSJacob Faibussowitsch   PetscCall(TSGetTimeError(ts, 0, &Yerr));
12369566063dSJacob Faibussowitsch   PetscCall(VecNorm(Yerr, NORM_2, &err_norm));
12379566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Yerr));
123863a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Estimated Error = %e.\n", (double)err_norm));
1239c4762a1bSJed Brown 
1240c4762a1bSJed Brown   /* Exact solution */
12419566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(Y, &Yex));
124263a3b9bcSJacob Faibussowitsch   if (PetscAbsReal(final_time - tfinal) > 2. * PETSC_MACHINE_EPSILON) {
12439566063dSJacob 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));
1244c4762a1bSJed Brown   }
12459566063dSJacob Faibussowitsch   PetscCall(ExactSolution(Yex, &ptype[0], final_time, exact_flag));
1246c4762a1bSJed Brown 
1247c4762a1bSJed Brown   /* Calculate Error */
12489566063dSJacob Faibussowitsch   PetscCall(VecAYPX(Yex, -1.0, Y));
12499566063dSJacob Faibussowitsch   PetscCall(VecNorm(Yex, NORM_2, error));
1250c4762a1bSJed Brown   *error = PetscSqrtReal(((*error) * (*error)) / N);
1251c4762a1bSJed Brown 
1252c4762a1bSJed Brown   /* Clean up and finalize */
12539566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Jac));
12549566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
12559566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Yex));
12569566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&Y));
1257c4762a1bSJed Brown 
1258*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
1259c4762a1bSJed Brown }
1260c4762a1bSJed Brown 
1261d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
1262d71ae5a4SJacob Faibussowitsch {
1263c4762a1bSJed Brown   char        ptype[256] = "hull1972a1"; /* Problem specification                                */
1264c4762a1bSJed Brown   PetscInt    n_refine   = 1;            /* Number of refinement levels for convergence analysis */
1265c4762a1bSJed Brown   PetscReal   refine_fac = 2.0;          /* Refinement factor for dt                             */
1266c4762a1bSJed Brown   PetscReal   dt_initial = 0.01;         /* Initial default value of dt                          */
1267c4762a1bSJed Brown   PetscReal   dt;
1268c4762a1bSJed Brown   PetscReal   tfinal  = 20.0;   /* Final time for the time-integration                  */
1269c4762a1bSJed Brown   PetscInt    maxiter = 100000; /* Maximum number of time-integration iterations        */
1270c4762a1bSJed Brown   PetscReal  *error;            /* Array to store the errors for convergence analysis   */
1271c4762a1bSJed Brown   PetscMPIInt size;             /* No of processors                                     */
1272c4762a1bSJed Brown   PetscBool   flag;             /* Flag denoting availability of exact solution         */
127311cc89d2SBarry Smith   PetscInt    r, N;
1274c4762a1bSJed Brown 
1275c4762a1bSJed Brown   /* Initialize program */
1276327415f7SBarry Smith   PetscFunctionBeginUser;
127711cc89d2SBarry Smith   PetscCall(GetSize(&ptype[0], &N));
12789566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
1279c4762a1bSJed Brown 
1280c4762a1bSJed Brown   /* Check if running with only 1 proc */
12819566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
12823c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
1283c4762a1bSJed Brown 
1284d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "ex31", NULL);
12859566063dSJacob Faibussowitsch   PetscCall(PetscOptionsString("-problem", "Problem specification", "<hull1972a1>", ptype, ptype, sizeof(ptype), NULL));
12869566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-refinement_levels", "Number of refinement levels for convergence analysis", "<1>", n_refine, &n_refine, NULL));
12879566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-refinement_factor", "Refinement factor for dt", "<2.0>", refine_fac, &refine_fac, NULL));
12889566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-dt", "Time step size (for convergence analysis, initial time step)", "<0.01>", dt_initial, &dt_initial, NULL));
12899566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-final_time", "Final time for the time-integration", "<20.0>", tfinal, &tfinal, NULL));
1290d0609cedSBarry Smith   PetscOptionsEnd();
1291c4762a1bSJed Brown 
12929566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(n_refine, &error));
1293c4762a1bSJed Brown   for (r = 0, dt = dt_initial; r < n_refine; r++) {
1294c4762a1bSJed Brown     error[r] = 0;
1295c4762a1bSJed Brown     if (r > 0) dt /= refine_fac;
1296c4762a1bSJed Brown 
129711cc89d2SBarry Smith     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, N));
12989566063dSJacob Faibussowitsch     PetscCall(SolveODE(&ptype[0], dt, tfinal, maxiter, &error[r], &flag));
1299c4762a1bSJed Brown     if (flag) {
1300c4762a1bSJed Brown       /* If exact solution available for the specified ODE */
1301c4762a1bSJed Brown       if (r > 0) {
1302c4762a1bSJed Brown         PetscReal conv_rate = (PetscLogReal(error[r]) - PetscLogReal(error[r - 1])) / (-PetscLogReal(refine_fac));
13039566063dSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error           = %E,\tConvergence rate = %f.\n", (double)error[r], (double)conv_rate));
1304c4762a1bSJed Brown       } else {
130563a3b9bcSJacob Faibussowitsch         PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Error           = %E.\n", (double)error[r]));
1306c4762a1bSJed Brown       }
1307c4762a1bSJed Brown     }
1308c4762a1bSJed Brown   }
13099566063dSJacob Faibussowitsch   PetscCall(PetscFree(error));
13109566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1311b122ec5aSJacob Faibussowitsch   return 0;
1312c4762a1bSJed Brown }
1313c4762a1bSJed Brown 
1314c4762a1bSJed Brown /*TEST
1315c4762a1bSJed Brown 
1316c4762a1bSJed Brown     test:
1317c4762a1bSJed Brown       suffix: 2
1318c4762a1bSJed Brown       args: -ts_type glee -final_time 5 -ts_adapt_type none
1319c4762a1bSJed Brown       timeoutfactor: 3
1320c4762a1bSJed Brown       requires: !single
1321c4762a1bSJed Brown 
1322c4762a1bSJed Brown     test:
1323c4762a1bSJed Brown       suffix: 3
1324c4762a1bSJed 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
1325c4762a1bSJed Brown       timeoutfactor: 3
1326c4762a1bSJed Brown       requires: !single
1327c4762a1bSJed Brown 
1328c4762a1bSJed Brown     test:
1329c4762a1bSJed Brown       suffix: 4
1330c4762a1bSJed 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
1331c4762a1bSJed Brown       timeoutfactor: 3
1332c4762a1bSJed Brown       requires: !single !__float128
1333c4762a1bSJed Brown 
1334c4762a1bSJed Brown TEST*/
1335