xref: /petsc/src/ts/tutorials/ex31.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Solves the ordinary differential equations (IVPs) using explicit and implicit time-integration methods.\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown   Concepts:   TS
6*c4762a1bSJed Brown   Useful command line parameters:
7*c4762a1bSJed Brown   -problem <hull1972a1>: choose which problem to solve (see references
8*c4762a1bSJed Brown                       for complete listing of problems).
9*c4762a1bSJed Brown   -ts_type <euler>: specify time-integrator
10*c4762a1bSJed Brown   -ts_adapt_type <basic>: specify time-step adapting (none,basic,advanced)
11*c4762a1bSJed Brown   -refinement_levels <1>: number of refinement levels for convergence analysis
12*c4762a1bSJed Brown   -refinement_factor <2.0>: factor to refine time step size by for convergence analysis
13*c4762a1bSJed Brown   -dt <0.01>: specify time step (initial time step for convergence analysis)
14*c4762a1bSJed Brown 
15*c4762a1bSJed Brown */
16*c4762a1bSJed Brown 
17*c4762a1bSJed Brown /*
18*c4762a1bSJed Brown List of cases and their names in the code:-
19*c4762a1bSJed Brown   From Hull, T.E., Enright, W.H., Fellen, B.M., and Sedgwick, A.E.,
20*c4762a1bSJed Brown       "Comparing Numerical Methods for Ordinary Differential
21*c4762a1bSJed Brown        Equations", SIAM J. Numer. Anal., 9(4), 1972, pp. 603 - 635
22*c4762a1bSJed Brown     A1 -> "hull1972a1" (exact solution available)
23*c4762a1bSJed Brown     A2 -> "hull1972a2" (exact solution available)
24*c4762a1bSJed Brown     A3 -> "hull1972a3" (exact solution available)
25*c4762a1bSJed Brown     A4 -> "hull1972a4" (exact solution available)
26*c4762a1bSJed Brown     A5 -> "hull1972a5"
27*c4762a1bSJed Brown     B1 -> "hull1972b1"
28*c4762a1bSJed Brown     B2 -> "hull1972b2"
29*c4762a1bSJed Brown     B3 -> "hull1972b3"
30*c4762a1bSJed Brown     B4 -> "hull1972b4"
31*c4762a1bSJed Brown     B5 -> "hull1972b5"
32*c4762a1bSJed Brown     C1 -> "hull1972c1"
33*c4762a1bSJed Brown     C2 -> "hull1972c2"
34*c4762a1bSJed Brown     C3 -> "hull1972c3"
35*c4762a1bSJed Brown     C4 -> "hull1972c4"
36*c4762a1bSJed Brown 
37*c4762a1bSJed Brown  From Constantinescu, E. "Estimating Global Errors in Time Stepping" ArXiv e-prints,
38*c4762a1bSJed Brown        https://arxiv.org/abs/1503.05166, 2016
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown     Kulikov2013I -> "kulik2013i"
41*c4762a1bSJed Brown 
42*c4762a1bSJed Brown */
43*c4762a1bSJed Brown 
44*c4762a1bSJed Brown #include <petscts.h>
45*c4762a1bSJed Brown 
46*c4762a1bSJed Brown /* Function declarations */
47*c4762a1bSJed Brown PetscErrorCode (*RHSFunction) (TS,PetscReal,Vec,Vec,void*);
48*c4762a1bSJed Brown PetscErrorCode (*RHSJacobian) (TS,PetscReal,Vec,Mat,Mat,void*);
49*c4762a1bSJed Brown PetscErrorCode (*IFunction)   (TS,PetscReal,Vec,Vec,Vec,void*);
50*c4762a1bSJed Brown PetscErrorCode (*IJacobian)   (TS,PetscReal,Vec,Vec,PetscReal,Mat,Mat,void*);
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown /* Returns the size of the system of equations depending on problem specification */
53*c4762a1bSJed Brown PetscInt GetSize(const char *p)
54*c4762a1bSJed Brown {
55*c4762a1bSJed Brown   PetscFunctionBegin;
56*c4762a1bSJed Brown   if      ((!strcmp(p,"hull1972a1"))
57*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972a2"))
58*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972a3"))
59*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972a4"))
60*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972a5")) )  PetscFunctionReturn(1);
61*c4762a1bSJed Brown   else if  (!strcmp(p,"hull1972b1")  )  PetscFunctionReturn(2);
62*c4762a1bSJed Brown   else if ((!strcmp(p,"hull1972b2"))
63*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972b3"))
64*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972b4"))
65*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972b5")) )  PetscFunctionReturn(3);
66*c4762a1bSJed Brown   else if ((!strcmp(p,"kulik2013i")) )  PetscFunctionReturn(4);
67*c4762a1bSJed Brown   else if ((!strcmp(p,"hull1972c1"))
68*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972c2"))
69*c4762a1bSJed Brown          ||(!strcmp(p,"hull1972c3")) )  PetscFunctionReturn(10);
70*c4762a1bSJed Brown   else if  (!strcmp(p,"hull1972c4")  )  PetscFunctionReturn(51);
71*c4762a1bSJed Brown   else                                  PetscFunctionReturn(-1);
72*c4762a1bSJed Brown }
73*c4762a1bSJed Brown 
74*c4762a1bSJed Brown /****************************************************************/
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown /* Problem specific functions */
77*c4762a1bSJed Brown 
78*c4762a1bSJed Brown /* Hull, 1972, Problem A1 */
79*c4762a1bSJed Brown 
80*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec F, void *s)
81*c4762a1bSJed Brown {
82*c4762a1bSJed Brown   PetscErrorCode    ierr;
83*c4762a1bSJed Brown   PetscScalar       *f;
84*c4762a1bSJed Brown   const PetscScalar *y;
85*c4762a1bSJed Brown 
86*c4762a1bSJed Brown   PetscFunctionBegin;
87*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
88*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
89*c4762a1bSJed Brown   f[0] = -y[0];
90*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
91*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
92*c4762a1bSJed Brown   PetscFunctionReturn(0);
93*c4762a1bSJed Brown }
94*c4762a1bSJed Brown 
95*c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
96*c4762a1bSJed Brown {
97*c4762a1bSJed Brown   PetscErrorCode    ierr;
98*c4762a1bSJed Brown   const PetscScalar *y;
99*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
100*c4762a1bSJed Brown   PetscScalar       value = -1.0;
101*c4762a1bSJed Brown 
102*c4762a1bSJed Brown   PetscFunctionBegin;
103*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
104*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
105*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
106*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
107*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
108*c4762a1bSJed Brown   PetscFunctionReturn(0);
109*c4762a1bSJed Brown }
110*c4762a1bSJed Brown 
111*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
112*c4762a1bSJed Brown {
113*c4762a1bSJed Brown   PetscErrorCode    ierr;
114*c4762a1bSJed Brown   const PetscScalar *y;
115*c4762a1bSJed Brown   PetscScalar       *f;
116*c4762a1bSJed Brown 
117*c4762a1bSJed Brown   PetscFunctionBegin;
118*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
119*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
120*c4762a1bSJed Brown   f[0] = -y[0];
121*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
122*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
123*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
124*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
125*c4762a1bSJed Brown   PetscFunctionReturn(0);
126*c4762a1bSJed Brown }
127*c4762a1bSJed Brown 
128*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
129*c4762a1bSJed Brown {
130*c4762a1bSJed Brown   PetscErrorCode    ierr;
131*c4762a1bSJed Brown   const PetscScalar *y;
132*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
133*c4762a1bSJed Brown   PetscScalar       value = a - 1.0;
134*c4762a1bSJed Brown 
135*c4762a1bSJed Brown   PetscFunctionBegin;
136*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
137*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
138*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
139*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
140*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
141*c4762a1bSJed Brown   PetscFunctionReturn(0);
142*c4762a1bSJed Brown }
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown /* Hull, 1972, Problem A2 */
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec F, void *s)
147*c4762a1bSJed Brown {
148*c4762a1bSJed Brown   PetscErrorCode    ierr;
149*c4762a1bSJed Brown   const PetscScalar *y;
150*c4762a1bSJed Brown   PetscScalar       *f;
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown   PetscFunctionBegin;
153*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
155*c4762a1bSJed Brown   f[0] = -0.5*y[0]*y[0]*y[0];
156*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
157*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
158*c4762a1bSJed Brown   PetscFunctionReturn(0);
159*c4762a1bSJed Brown }
160*c4762a1bSJed Brown 
161*c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
162*c4762a1bSJed Brown {
163*c4762a1bSJed Brown   PetscErrorCode    ierr;
164*c4762a1bSJed Brown   const PetscScalar *y;
165*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
166*c4762a1bSJed Brown   PetscScalar       value;
167*c4762a1bSJed Brown 
168*c4762a1bSJed Brown   PetscFunctionBegin;
169*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
170*c4762a1bSJed Brown   value = -0.5*3.0*y[0]*y[0];
171*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
174*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
175*c4762a1bSJed Brown   PetscFunctionReturn(0);
176*c4762a1bSJed Brown }
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
179*c4762a1bSJed Brown {
180*c4762a1bSJed Brown   PetscErrorCode    ierr;
181*c4762a1bSJed Brown   PetscScalar       *f;
182*c4762a1bSJed Brown   const PetscScalar *y;
183*c4762a1bSJed Brown 
184*c4762a1bSJed Brown   PetscFunctionBegin;
185*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
186*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
187*c4762a1bSJed Brown   f[0] = -0.5*y[0]*y[0]*y[0];
188*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
189*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
190*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
191*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
192*c4762a1bSJed Brown   PetscFunctionReturn(0);
193*c4762a1bSJed Brown }
194*c4762a1bSJed Brown 
195*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
196*c4762a1bSJed Brown {
197*c4762a1bSJed Brown   PetscErrorCode    ierr;
198*c4762a1bSJed Brown   const PetscScalar *y;
199*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
200*c4762a1bSJed Brown   PetscScalar       value;
201*c4762a1bSJed Brown 
202*c4762a1bSJed Brown   PetscFunctionBegin;
203*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
204*c4762a1bSJed Brown   value = a + 0.5*3.0*y[0]*y[0];
205*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
206*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
207*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
208*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
209*c4762a1bSJed Brown   PetscFunctionReturn(0);
210*c4762a1bSJed Brown }
211*c4762a1bSJed Brown 
212*c4762a1bSJed Brown /* Hull, 1972, Problem A3 */
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec F, void *s)
215*c4762a1bSJed Brown {
216*c4762a1bSJed Brown   PetscErrorCode    ierr;
217*c4762a1bSJed Brown   const PetscScalar *y;
218*c4762a1bSJed Brown   PetscScalar       *f;
219*c4762a1bSJed Brown 
220*c4762a1bSJed Brown   PetscFunctionBegin;
221*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
222*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
223*c4762a1bSJed Brown   f[0] = y[0]*PetscCosReal(t);
224*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
225*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
226*c4762a1bSJed Brown   PetscFunctionReturn(0);
227*c4762a1bSJed Brown }
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
230*c4762a1bSJed Brown {
231*c4762a1bSJed Brown   PetscErrorCode    ierr;
232*c4762a1bSJed Brown   const PetscScalar *y;
233*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
234*c4762a1bSJed Brown   PetscScalar       value = PetscCosReal(t);
235*c4762a1bSJed Brown 
236*c4762a1bSJed Brown   PetscFunctionBegin;
237*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
238*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
239*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
240*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
241*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
242*c4762a1bSJed Brown   PetscFunctionReturn(0);
243*c4762a1bSJed Brown }
244*c4762a1bSJed Brown 
245*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
246*c4762a1bSJed Brown {
247*c4762a1bSJed Brown   PetscErrorCode    ierr;
248*c4762a1bSJed Brown   PetscScalar       *f;
249*c4762a1bSJed Brown   const PetscScalar *y;
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown   PetscFunctionBegin;
252*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
253*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
254*c4762a1bSJed Brown   f[0] = y[0]*PetscCosReal(t);
255*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
256*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
257*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
258*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
259*c4762a1bSJed Brown   PetscFunctionReturn(0);
260*c4762a1bSJed Brown }
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
263*c4762a1bSJed Brown {
264*c4762a1bSJed Brown   PetscErrorCode    ierr;
265*c4762a1bSJed Brown   const PetscScalar *y;
266*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
267*c4762a1bSJed Brown   PetscScalar       value = a - PetscCosReal(t);
268*c4762a1bSJed Brown 
269*c4762a1bSJed Brown   PetscFunctionBegin;
270*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
271*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
272*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
273*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
274*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
275*c4762a1bSJed Brown   PetscFunctionReturn(0);
276*c4762a1bSJed Brown }
277*c4762a1bSJed Brown 
278*c4762a1bSJed Brown /* Hull, 1972, Problem A4 */
279*c4762a1bSJed Brown 
280*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec F, void *s)
281*c4762a1bSJed Brown {
282*c4762a1bSJed Brown   PetscErrorCode    ierr;
283*c4762a1bSJed Brown   PetscScalar       *f;
284*c4762a1bSJed Brown   const PetscScalar *y;
285*c4762a1bSJed Brown 
286*c4762a1bSJed Brown   PetscFunctionBegin;
287*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
288*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
289*c4762a1bSJed Brown   f[0] = (0.25*y[0])*(1.0-0.05*y[0]);
290*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
291*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
292*c4762a1bSJed Brown   PetscFunctionReturn(0);
293*c4762a1bSJed Brown }
294*c4762a1bSJed Brown 
295*c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
296*c4762a1bSJed Brown {
297*c4762a1bSJed Brown   PetscErrorCode    ierr;
298*c4762a1bSJed Brown   const PetscScalar *y;
299*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
300*c4762a1bSJed Brown   PetscScalar       value;
301*c4762a1bSJed Brown 
302*c4762a1bSJed Brown   PetscFunctionBegin;
303*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
304*c4762a1bSJed Brown   value = 0.25*(1.0-0.05*y[0]) - (0.25*y[0])*0.05;
305*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
306*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
307*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
308*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
309*c4762a1bSJed Brown   PetscFunctionReturn(0);
310*c4762a1bSJed Brown }
311*c4762a1bSJed Brown 
312*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
313*c4762a1bSJed Brown {
314*c4762a1bSJed Brown   PetscErrorCode    ierr;
315*c4762a1bSJed Brown   PetscScalar       *f;
316*c4762a1bSJed Brown   const PetscScalar *y;
317*c4762a1bSJed Brown 
318*c4762a1bSJed Brown   PetscFunctionBegin;
319*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
320*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
321*c4762a1bSJed Brown   f[0] = (0.25*y[0])*(1.0-0.05*y[0]);
322*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
323*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
324*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
325*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
326*c4762a1bSJed Brown   PetscFunctionReturn(0);
327*c4762a1bSJed Brown }
328*c4762a1bSJed Brown 
329*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
330*c4762a1bSJed Brown {
331*c4762a1bSJed Brown   PetscErrorCode    ierr;
332*c4762a1bSJed Brown   const PetscScalar *y;
333*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
334*c4762a1bSJed Brown   PetscScalar       value;
335*c4762a1bSJed Brown 
336*c4762a1bSJed Brown   PetscFunctionBegin;
337*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
338*c4762a1bSJed Brown   value = a - 0.25*(1.0-0.05*y[0]) + (0.25*y[0])*0.05;
339*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
341*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
342*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
343*c4762a1bSJed Brown   PetscFunctionReturn(0);
344*c4762a1bSJed Brown }
345*c4762a1bSJed Brown 
346*c4762a1bSJed Brown /* Hull, 1972, Problem A5 */
347*c4762a1bSJed Brown 
348*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec F, void *s)
349*c4762a1bSJed Brown {
350*c4762a1bSJed Brown   PetscErrorCode    ierr;
351*c4762a1bSJed Brown   PetscScalar       *f;
352*c4762a1bSJed Brown   const PetscScalar *y;
353*c4762a1bSJed Brown 
354*c4762a1bSJed Brown   PetscFunctionBegin;
355*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
356*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
357*c4762a1bSJed Brown   f[0] = (y[0]-t)/(y[0]+t);
358*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
359*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
360*c4762a1bSJed Brown   PetscFunctionReturn(0);
361*c4762a1bSJed Brown }
362*c4762a1bSJed Brown 
363*c4762a1bSJed Brown PetscErrorCode RHSJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
364*c4762a1bSJed Brown {
365*c4762a1bSJed Brown   PetscErrorCode    ierr;
366*c4762a1bSJed Brown   const PetscScalar *y;
367*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
368*c4762a1bSJed Brown   PetscScalar       value;
369*c4762a1bSJed Brown 
370*c4762a1bSJed Brown   PetscFunctionBegin;
371*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
372*c4762a1bSJed Brown   value = 2*t/((t+y[0])*(t+y[0]));
373*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
374*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
375*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
376*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
377*c4762a1bSJed Brown   PetscFunctionReturn(0);
378*c4762a1bSJed Brown }
379*c4762a1bSJed Brown 
380*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
381*c4762a1bSJed Brown {
382*c4762a1bSJed Brown   PetscErrorCode    ierr;
383*c4762a1bSJed Brown   PetscScalar       *f;
384*c4762a1bSJed Brown   const PetscScalar *y;
385*c4762a1bSJed Brown 
386*c4762a1bSJed Brown   PetscFunctionBegin;
387*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
388*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
389*c4762a1bSJed Brown   f[0] = (y[0]-t)/(y[0]+t);
390*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
391*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
392*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
393*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
394*c4762a1bSJed Brown   PetscFunctionReturn(0);
395*c4762a1bSJed Brown }
396*c4762a1bSJed Brown 
397*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972A5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
398*c4762a1bSJed Brown {
399*c4762a1bSJed Brown   PetscErrorCode    ierr;
400*c4762a1bSJed Brown   const PetscScalar *y;
401*c4762a1bSJed Brown   PetscInt          row = 0,col = 0;
402*c4762a1bSJed Brown   PetscScalar       value;
403*c4762a1bSJed Brown 
404*c4762a1bSJed Brown   PetscFunctionBegin;
405*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
406*c4762a1bSJed Brown   value = a - 2*t/((t+y[0])*(t+y[0]));
407*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&row,1,&col,&value,INSERT_VALUES);CHKERRQ(ierr);
408*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
409*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
410*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
411*c4762a1bSJed Brown   PetscFunctionReturn(0);
412*c4762a1bSJed Brown }
413*c4762a1bSJed Brown 
414*c4762a1bSJed Brown /* Hull, 1972, Problem B1 */
415*c4762a1bSJed Brown 
416*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec F, void *s)
417*c4762a1bSJed Brown {
418*c4762a1bSJed Brown   PetscErrorCode    ierr;
419*c4762a1bSJed Brown   PetscScalar       *f;
420*c4762a1bSJed Brown   const PetscScalar *y;
421*c4762a1bSJed Brown 
422*c4762a1bSJed Brown   PetscFunctionBegin;
423*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
424*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
425*c4762a1bSJed Brown   f[0] = 2.0*(y[0] - y[0]*y[1]);
426*c4762a1bSJed Brown   f[1] = -(y[1]-y[0]*y[1]);
427*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
428*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
429*c4762a1bSJed Brown   PetscFunctionReturn(0);
430*c4762a1bSJed Brown }
431*c4762a1bSJed Brown 
432*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
433*c4762a1bSJed Brown {
434*c4762a1bSJed Brown   PetscErrorCode    ierr;
435*c4762a1bSJed Brown   PetscScalar       *f;
436*c4762a1bSJed Brown   const PetscScalar *y;
437*c4762a1bSJed Brown 
438*c4762a1bSJed Brown   PetscFunctionBegin;
439*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
440*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
441*c4762a1bSJed Brown   f[0] = 2.0*(y[0] - y[0]*y[1]);
442*c4762a1bSJed Brown   f[1] = -(y[1]-y[0]*y[1]);
443*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
444*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
445*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
446*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
447*c4762a1bSJed Brown   PetscFunctionReturn(0);
448*c4762a1bSJed Brown }
449*c4762a1bSJed Brown 
450*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
451*c4762a1bSJed Brown {
452*c4762a1bSJed Brown   PetscErrorCode    ierr;
453*c4762a1bSJed Brown   const PetscScalar *y;
454*c4762a1bSJed Brown   PetscInt          row[2] = {0,1};
455*c4762a1bSJed Brown   PetscScalar       value[2][2];
456*c4762a1bSJed Brown 
457*c4762a1bSJed Brown   PetscFunctionBegin;
458*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
459*c4762a1bSJed Brown   value[0][0] = a - 2.0*(1.0-y[1]);    value[0][1] = 2.0*y[0];
460*c4762a1bSJed Brown   value[1][0] = -y[1];                 value[1][1] = a + 1.0 - y[0];
461*c4762a1bSJed Brown   ierr = MatSetValues(A,2,&row[0],2,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr);
462*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
463*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
464*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
465*c4762a1bSJed Brown   PetscFunctionReturn(0);
466*c4762a1bSJed Brown }
467*c4762a1bSJed Brown 
468*c4762a1bSJed Brown /* Hull, 1972, Problem B2 */
469*c4762a1bSJed Brown 
470*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec F, void *s)
471*c4762a1bSJed Brown {
472*c4762a1bSJed Brown   PetscErrorCode    ierr;
473*c4762a1bSJed Brown   PetscScalar       *f;
474*c4762a1bSJed Brown   const PetscScalar *y;
475*c4762a1bSJed Brown 
476*c4762a1bSJed Brown   PetscFunctionBegin;
477*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
478*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
479*c4762a1bSJed Brown   f[0] = -y[0] + y[1];
480*c4762a1bSJed Brown   f[1] = y[0] - 2.0*y[1] + y[2];
481*c4762a1bSJed Brown   f[2] = y[1] - y[2];
482*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
483*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
484*c4762a1bSJed Brown   PetscFunctionReturn(0);
485*c4762a1bSJed Brown }
486*c4762a1bSJed Brown 
487*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
488*c4762a1bSJed Brown {
489*c4762a1bSJed Brown   PetscErrorCode    ierr;
490*c4762a1bSJed Brown   PetscScalar       *f;
491*c4762a1bSJed Brown   const PetscScalar *y;
492*c4762a1bSJed Brown 
493*c4762a1bSJed Brown   PetscFunctionBegin;
494*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
495*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
496*c4762a1bSJed Brown   f[0] = -y[0] + y[1];
497*c4762a1bSJed Brown   f[1] = y[0] - 2.0*y[1] + y[2];
498*c4762a1bSJed Brown   f[2] = y[1] - y[2];
499*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
500*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
501*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
502*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
503*c4762a1bSJed Brown   PetscFunctionReturn(0);
504*c4762a1bSJed Brown }
505*c4762a1bSJed Brown 
506*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
507*c4762a1bSJed Brown {
508*c4762a1bSJed Brown   PetscErrorCode    ierr;
509*c4762a1bSJed Brown   const PetscScalar *y;
510*c4762a1bSJed Brown   PetscInt          row[3] = {0,1,2};
511*c4762a1bSJed Brown   PetscScalar       value[3][3];
512*c4762a1bSJed Brown 
513*c4762a1bSJed Brown   PetscFunctionBegin;
514*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
515*c4762a1bSJed Brown   value[0][0] = a + 1.0;  value[0][1] = -1.0;     value[0][2] = 0;
516*c4762a1bSJed Brown   value[1][0] = -1.0;     value[1][1] = a + 2.0;  value[1][2] = -1.0;
517*c4762a1bSJed Brown   value[2][0] = 0;        value[2][1] = -1.0;     value[2][2] = a + 1.0;
518*c4762a1bSJed Brown   ierr = MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr);
519*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
520*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
521*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
522*c4762a1bSJed Brown   PetscFunctionReturn(0);
523*c4762a1bSJed Brown }
524*c4762a1bSJed Brown 
525*c4762a1bSJed Brown /* Hull, 1972, Problem B3 */
526*c4762a1bSJed Brown 
527*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec F, void *s)
528*c4762a1bSJed Brown {
529*c4762a1bSJed Brown   PetscErrorCode    ierr;
530*c4762a1bSJed Brown   PetscScalar       *f;
531*c4762a1bSJed Brown   const PetscScalar *y;
532*c4762a1bSJed Brown 
533*c4762a1bSJed Brown   PetscFunctionBegin;
534*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
535*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
536*c4762a1bSJed Brown   f[0] = -y[0];
537*c4762a1bSJed Brown   f[1] = y[0] - y[1]*y[1];
538*c4762a1bSJed Brown   f[2] = y[1]*y[1];
539*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
540*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
541*c4762a1bSJed Brown   PetscFunctionReturn(0);
542*c4762a1bSJed Brown }
543*c4762a1bSJed Brown 
544*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
545*c4762a1bSJed Brown {
546*c4762a1bSJed Brown   PetscErrorCode    ierr;
547*c4762a1bSJed Brown   PetscScalar       *f;
548*c4762a1bSJed Brown   const PetscScalar *y;
549*c4762a1bSJed Brown 
550*c4762a1bSJed Brown   PetscFunctionBegin;
551*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
552*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
553*c4762a1bSJed Brown   f[0] = -y[0];
554*c4762a1bSJed Brown   f[1] = y[0] - y[1]*y[1];
555*c4762a1bSJed Brown   f[2] = y[1]*y[1];
556*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
557*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
558*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
559*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
560*c4762a1bSJed Brown   PetscFunctionReturn(0);
561*c4762a1bSJed Brown }
562*c4762a1bSJed Brown 
563*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B3(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
564*c4762a1bSJed Brown {
565*c4762a1bSJed Brown   PetscErrorCode    ierr;
566*c4762a1bSJed Brown   const PetscScalar *y;
567*c4762a1bSJed Brown   PetscInt          row[3] = {0,1,2};
568*c4762a1bSJed Brown   PetscScalar       value[3][3];
569*c4762a1bSJed Brown 
570*c4762a1bSJed Brown   PetscFunctionBegin;
571*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
572*c4762a1bSJed Brown   value[0][0] = a + 1.0; value[0][1] = 0;             value[0][2] = 0;
573*c4762a1bSJed Brown   value[1][0] = -1.0;    value[1][1] = a + 2.0*y[1];  value[1][2] = 0;
574*c4762a1bSJed Brown   value[2][0] = 0;       value[2][1] = -2.0*y[1];     value[2][2] = a;
575*c4762a1bSJed Brown   ierr = MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr);
576*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
577*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
578*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
579*c4762a1bSJed Brown   PetscFunctionReturn(0);
580*c4762a1bSJed Brown }
581*c4762a1bSJed Brown 
582*c4762a1bSJed Brown /* Hull, 1972, Problem B4 */
583*c4762a1bSJed Brown 
584*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec F, void *s)
585*c4762a1bSJed Brown {
586*c4762a1bSJed Brown   PetscErrorCode    ierr;
587*c4762a1bSJed Brown   PetscScalar       *f;
588*c4762a1bSJed Brown   const PetscScalar *y;
589*c4762a1bSJed Brown 
590*c4762a1bSJed Brown   PetscFunctionBegin;
591*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
592*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
593*c4762a1bSJed Brown   f[0] = -y[1] - y[0]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]);
594*c4762a1bSJed Brown   f[1] =  y[0] - y[1]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]);
595*c4762a1bSJed Brown   f[2] = y[0]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]);
596*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
597*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
598*c4762a1bSJed Brown   PetscFunctionReturn(0);
599*c4762a1bSJed Brown }
600*c4762a1bSJed Brown 
601*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
602*c4762a1bSJed Brown {
603*c4762a1bSJed Brown   PetscErrorCode    ierr;
604*c4762a1bSJed Brown   PetscScalar       *f;
605*c4762a1bSJed Brown   const PetscScalar *y;
606*c4762a1bSJed Brown 
607*c4762a1bSJed Brown   PetscFunctionBegin;
608*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
609*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
610*c4762a1bSJed Brown   f[0] = -y[1] - y[0]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]);
611*c4762a1bSJed Brown   f[1] =  y[0] - y[1]*y[2]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]);
612*c4762a1bSJed Brown   f[2] = y[0]/PetscSqrtScalar(y[0]*y[0]+y[1]*y[1]);
613*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
614*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
615*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
616*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
617*c4762a1bSJed Brown   PetscFunctionReturn(0);
618*c4762a1bSJed Brown }
619*c4762a1bSJed Brown 
620*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B4(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
621*c4762a1bSJed Brown {
622*c4762a1bSJed Brown   PetscErrorCode    ierr;
623*c4762a1bSJed Brown   const PetscScalar *y;
624*c4762a1bSJed Brown   PetscInt          row[3] = {0,1,2};
625*c4762a1bSJed Brown   PetscScalar       value[3][3],fac,fac2;
626*c4762a1bSJed Brown 
627*c4762a1bSJed Brown   PetscFunctionBegin;
628*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
629*c4762a1bSJed Brown   fac  = PetscPowScalar(y[0]*y[0]+y[1]*y[1],-1.5);
630*c4762a1bSJed Brown   fac2 = PetscPowScalar(y[0]*y[0]+y[1]*y[1],-0.5);
631*c4762a1bSJed Brown   value[0][0] = a + (y[1]*y[1]*y[2])*fac;
632*c4762a1bSJed Brown   value[0][1] = 1.0 - (y[0]*y[1]*y[2])*fac;
633*c4762a1bSJed Brown   value[0][2] = y[0]*fac2;
634*c4762a1bSJed Brown   value[1][0] = -1.0 - y[0]*y[1]*y[2]*fac;
635*c4762a1bSJed Brown   value[1][1] = a + y[0]*y[0]*y[2]*fac;
636*c4762a1bSJed Brown   value[1][2] = y[1]*fac2;
637*c4762a1bSJed Brown   value[2][0] = -y[1]*y[1]*fac;
638*c4762a1bSJed Brown   value[2][1] = y[0]*y[1]*fac;
639*c4762a1bSJed Brown   value[2][2] = a;
640*c4762a1bSJed Brown   ierr = MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr);
641*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
642*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
643*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
644*c4762a1bSJed Brown   PetscFunctionReturn(0);
645*c4762a1bSJed Brown }
646*c4762a1bSJed Brown 
647*c4762a1bSJed Brown /* Hull, 1972, Problem B5 */
648*c4762a1bSJed Brown 
649*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec F, void *s)
650*c4762a1bSJed Brown {
651*c4762a1bSJed Brown   PetscErrorCode    ierr;
652*c4762a1bSJed Brown   PetscScalar       *f;
653*c4762a1bSJed Brown   const PetscScalar *y;
654*c4762a1bSJed Brown 
655*c4762a1bSJed Brown   PetscFunctionBegin;
656*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
657*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
658*c4762a1bSJed Brown   f[0] = y[1]*y[2];
659*c4762a1bSJed Brown   f[1] = -y[0]*y[2];
660*c4762a1bSJed Brown   f[2] = -0.51*y[0]*y[1];
661*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
662*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
663*c4762a1bSJed Brown   PetscFunctionReturn(0);
664*c4762a1bSJed Brown }
665*c4762a1bSJed Brown 
666*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
667*c4762a1bSJed Brown {
668*c4762a1bSJed Brown   PetscErrorCode    ierr;
669*c4762a1bSJed Brown   PetscScalar       *f;
670*c4762a1bSJed Brown   const PetscScalar *y;
671*c4762a1bSJed Brown 
672*c4762a1bSJed Brown   PetscFunctionBegin;
673*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
674*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
675*c4762a1bSJed Brown   f[0] = y[1]*y[2];
676*c4762a1bSJed Brown   f[1] = -y[0]*y[2];
677*c4762a1bSJed Brown   f[2] = -0.51*y[0]*y[1];
678*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
679*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
680*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
681*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
682*c4762a1bSJed Brown   PetscFunctionReturn(0);
683*c4762a1bSJed Brown }
684*c4762a1bSJed Brown 
685*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972B5(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
686*c4762a1bSJed Brown {
687*c4762a1bSJed Brown   PetscErrorCode    ierr;
688*c4762a1bSJed Brown   const PetscScalar *y;
689*c4762a1bSJed Brown   PetscInt          row[3] = {0,1,2};
690*c4762a1bSJed Brown   PetscScalar       value[3][3];
691*c4762a1bSJed Brown 
692*c4762a1bSJed Brown   PetscFunctionBegin;
693*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
694*c4762a1bSJed Brown   value[0][0] = a;          value[0][1] = -y[2];      value[0][2] = -y[1];
695*c4762a1bSJed Brown   value[1][0] = y[2];       value[1][1] = a;          value[1][2] = y[0];
696*c4762a1bSJed Brown   value[2][0] = 0.51*y[1];  value[2][1] = 0.51*y[0];  value[2][2] = a;
697*c4762a1bSJed Brown   ierr = MatSetValues(A,3,&row[0],3,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr);
698*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
699*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
700*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
701*c4762a1bSJed Brown   PetscFunctionReturn(0);
702*c4762a1bSJed Brown }
703*c4762a1bSJed Brown 
704*c4762a1bSJed Brown 
705*c4762a1bSJed Brown /* Kulikov, 2013, Problem I */
706*c4762a1bSJed Brown 
707*c4762a1bSJed Brown PetscErrorCode RHSFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec F, void *s)
708*c4762a1bSJed Brown {
709*c4762a1bSJed Brown   PetscErrorCode    ierr;
710*c4762a1bSJed Brown   PetscScalar       *f;
711*c4762a1bSJed Brown   const PetscScalar *y;
712*c4762a1bSJed Brown 
713*c4762a1bSJed Brown   PetscFunctionBegin;
714*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
715*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
716*c4762a1bSJed Brown   f[0] = 2.*t*PetscPowScalar(y[1],1./5.)*y[3];
717*c4762a1bSJed Brown   f[1] = 10.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.));
718*c4762a1bSJed Brown   f[2] = 2.*t*y[3];
719*c4762a1bSJed Brown   f[3] = -2.*t*PetscLogScalar(y[0]);
720*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
721*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
722*c4762a1bSJed Brown   PetscFunctionReturn(0);
723*c4762a1bSJed Brown }
724*c4762a1bSJed Brown 
725*c4762a1bSJed Brown PetscErrorCode RHSJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Mat A, Mat B, void *s)
726*c4762a1bSJed Brown {
727*c4762a1bSJed Brown   PetscErrorCode    ierr;
728*c4762a1bSJed Brown   const PetscScalar *y;
729*c4762a1bSJed Brown   PetscInt          row[4] = {0,1,2,3};
730*c4762a1bSJed Brown   PetscScalar       value[4][4];
731*c4762a1bSJed Brown   PetscScalar       m1,m2;
732*c4762a1bSJed Brown   PetscFunctionBegin;
733*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
734*c4762a1bSJed Brown   m1=(2.*t*y[3])/(5.*PetscPowScalar(y[1],4./5.));
735*c4762a1bSJed Brown   m2=2.*t*PetscPowScalar(y[1],1./5.);
736*c4762a1bSJed Brown   value[0][0] = 0. ;        value[0][1] = m1; value[0][2] = 0.;  value[0][3] = m2;
737*c4762a1bSJed Brown   m1=50.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.));
738*c4762a1bSJed Brown   m2=10.*t*PetscExpScalar(5.0*(y[2]-1.));
739*c4762a1bSJed Brown   value[1][0] = 0.;        value[1][1] = 0. ; value[1][2] = m1; value[1][3] = m2;
740*c4762a1bSJed Brown   value[2][0] = 0.;        value[2][1] = 0.;  value[2][2] = 0.; value[2][3] = 2*t;
741*c4762a1bSJed Brown   value[3][0] = -2.*t/y[0];value[3][1] = 0.;  value[3][2] = 0.; value[3][3] = 0.;
742*c4762a1bSJed Brown   ierr = MatSetValues(A,4,&row[0],4,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr);
743*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
744*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
745*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
746*c4762a1bSJed Brown   PetscFunctionReturn(0);
747*c4762a1bSJed Brown }
748*c4762a1bSJed Brown 
749*c4762a1bSJed Brown PetscErrorCode IFunction_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
750*c4762a1bSJed Brown {
751*c4762a1bSJed Brown   PetscErrorCode    ierr;
752*c4762a1bSJed Brown   PetscScalar       *f;
753*c4762a1bSJed Brown   const PetscScalar *y;
754*c4762a1bSJed Brown 
755*c4762a1bSJed Brown   PetscFunctionBegin;
756*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
757*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
758*c4762a1bSJed Brown   f[0] = 2.*t*PetscPowScalar(y[1],1./5.)*y[3];
759*c4762a1bSJed Brown   f[1] = 10.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.));
760*c4762a1bSJed Brown   f[2] = 2.*t*y[3];
761*c4762a1bSJed Brown   f[3] = -2.*t*PetscLogScalar(y[0]);
762*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
763*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
764*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
765*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
766*c4762a1bSJed Brown   PetscFunctionReturn(0);
767*c4762a1bSJed Brown }
768*c4762a1bSJed Brown 
769*c4762a1bSJed Brown PetscErrorCode IJacobian_Kulikov2013I(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
770*c4762a1bSJed Brown {
771*c4762a1bSJed Brown   PetscErrorCode    ierr;
772*c4762a1bSJed Brown   const PetscScalar *y;
773*c4762a1bSJed Brown   PetscInt          row[4] = {0,1,2,3};
774*c4762a1bSJed Brown   PetscScalar       value[4][4];
775*c4762a1bSJed Brown   PetscScalar       m1,m2;
776*c4762a1bSJed Brown 
777*c4762a1bSJed Brown   PetscFunctionBegin;
778*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
779*c4762a1bSJed Brown   m1=(2.*t*y[3])/(5.*PetscPowScalar(y[1],4./5.));
780*c4762a1bSJed Brown   m2=2.*t*PetscPowScalar(y[1],1./5.);
781*c4762a1bSJed Brown   value[0][0] = a ;        value[0][1] = m1;  value[0][2] = 0.; value[0][3] = m2;
782*c4762a1bSJed Brown   m1=50.*t*y[3]*PetscExpScalar(5.0*(y[2]-1.));
783*c4762a1bSJed Brown   m2=10.*t*PetscExpScalar(5.0*(y[2]-1.));
784*c4762a1bSJed Brown   value[1][0] = 0.;        value[1][1] = a ;  value[1][2] = m1; value[1][3] = m2;
785*c4762a1bSJed Brown   value[2][0] = 0.;        value[2][1] = 0.;  value[2][2] = a;  value[2][3] = 2*t;
786*c4762a1bSJed Brown   value[3][0] = -2.*t/y[0];value[3][1] = 0.;  value[3][2] = 0.; value[3][3] = a;
787*c4762a1bSJed Brown   ierr = MatSetValues(A,4,&row[0],4,&row[0],&value[0][0],INSERT_VALUES);CHKERRQ(ierr);
788*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
789*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
790*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
791*c4762a1bSJed Brown   PetscFunctionReturn(0);
792*c4762a1bSJed Brown }
793*c4762a1bSJed Brown 
794*c4762a1bSJed Brown 
795*c4762a1bSJed Brown /* Hull, 1972, Problem C1 */
796*c4762a1bSJed Brown 
797*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec F, void *s)
798*c4762a1bSJed Brown {
799*c4762a1bSJed Brown   PetscErrorCode    ierr;
800*c4762a1bSJed Brown   PetscScalar       *f;
801*c4762a1bSJed Brown   const PetscScalar *y;
802*c4762a1bSJed Brown   PetscInt          N,i;
803*c4762a1bSJed Brown 
804*c4762a1bSJed Brown   PetscFunctionBegin;
805*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
806*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
807*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
808*c4762a1bSJed Brown   f[0] = -y[0];
809*c4762a1bSJed Brown   for (i = 1; i < N-1; i++) {
810*c4762a1bSJed Brown     f[i] = y[i-1] - y[i];
811*c4762a1bSJed Brown   }
812*c4762a1bSJed Brown   f[N-1] = y[N-2];
813*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
814*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
815*c4762a1bSJed Brown   PetscFunctionReturn(0);
816*c4762a1bSJed Brown }
817*c4762a1bSJed Brown 
818*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
819*c4762a1bSJed Brown {
820*c4762a1bSJed Brown   PetscErrorCode    ierr;
821*c4762a1bSJed Brown   PetscScalar       *f;
822*c4762a1bSJed Brown   const PetscScalar *y;
823*c4762a1bSJed Brown   PetscInt          N,i;
824*c4762a1bSJed Brown 
825*c4762a1bSJed Brown   PetscFunctionBegin;
826*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
827*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
828*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
829*c4762a1bSJed Brown   f[0] = -y[0];
830*c4762a1bSJed Brown   for (i = 1; i < N-1; i++) {
831*c4762a1bSJed Brown     f[i] = y[i-1] - y[i];
832*c4762a1bSJed Brown   }
833*c4762a1bSJed Brown   f[N-1] = y[N-2];
834*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
835*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
836*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
837*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
838*c4762a1bSJed Brown   PetscFunctionReturn(0);
839*c4762a1bSJed Brown }
840*c4762a1bSJed Brown 
841*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C1(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
842*c4762a1bSJed Brown {
843*c4762a1bSJed Brown   PetscErrorCode    ierr;
844*c4762a1bSJed Brown   const PetscScalar *y;
845*c4762a1bSJed Brown   PetscInt          N,i,col[2];
846*c4762a1bSJed Brown   PetscScalar       value[2];
847*c4762a1bSJed Brown 
848*c4762a1bSJed Brown   PetscFunctionBegin;
849*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
850*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
851*c4762a1bSJed Brown   i = 0;
852*c4762a1bSJed Brown   value[0] = a+1; col[0] = 0;
853*c4762a1bSJed Brown   value[1] =  0;  col[1] = 1;
854*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
855*c4762a1bSJed Brown   for (i = 0; i < N; i++) {
856*c4762a1bSJed Brown     value[0] =  -1; col[0] = i-1;
857*c4762a1bSJed Brown     value[1] = a+1; col[1] = i;
858*c4762a1bSJed Brown     ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
859*c4762a1bSJed Brown   }
860*c4762a1bSJed Brown   i = N-1;
861*c4762a1bSJed Brown   value[0] = -1;  col[0] = N-2;
862*c4762a1bSJed Brown   value[1] = a;   col[1] = N-1;
863*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
864*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
865*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
866*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
867*c4762a1bSJed Brown   PetscFunctionReturn(0);
868*c4762a1bSJed Brown }
869*c4762a1bSJed Brown 
870*c4762a1bSJed Brown /* Hull, 1972, Problem C2 */
871*c4762a1bSJed Brown 
872*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec F, void *s)
873*c4762a1bSJed Brown {
874*c4762a1bSJed Brown   PetscErrorCode    ierr;
875*c4762a1bSJed Brown   const PetscScalar *y;
876*c4762a1bSJed Brown   PetscScalar       *f;
877*c4762a1bSJed Brown   PetscInt          N,i;
878*c4762a1bSJed Brown 
879*c4762a1bSJed Brown   PetscFunctionBegin;
880*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
881*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
882*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
883*c4762a1bSJed Brown   f[0] = -y[0];
884*c4762a1bSJed Brown   for (i = 1; i < N-1; i++) {
885*c4762a1bSJed Brown     f[i] = (PetscReal)i*y[i-1] - (PetscReal)(i+1)*y[i];
886*c4762a1bSJed Brown   }
887*c4762a1bSJed Brown   f[N-1] = (PetscReal)(N-1)*y[N-2];
888*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
889*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
890*c4762a1bSJed Brown   PetscFunctionReturn(0);
891*c4762a1bSJed Brown }
892*c4762a1bSJed Brown 
893*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
894*c4762a1bSJed Brown {
895*c4762a1bSJed Brown   PetscErrorCode    ierr;
896*c4762a1bSJed Brown   PetscScalar       *f;
897*c4762a1bSJed Brown   const PetscScalar *y;
898*c4762a1bSJed Brown   PetscInt          N,i;
899*c4762a1bSJed Brown 
900*c4762a1bSJed Brown   PetscFunctionBegin;
901*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
902*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
903*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
904*c4762a1bSJed Brown   f[0] = -y[0];
905*c4762a1bSJed Brown   for (i = 1; i < N-1; i++) {
906*c4762a1bSJed Brown     f[i] = (PetscReal)i*y[i-1] - (PetscReal)(i+1)*y[i];
907*c4762a1bSJed Brown   }
908*c4762a1bSJed Brown   f[N-1] = (PetscReal)(N-1)*y[N-2];
909*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
910*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
911*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
912*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
913*c4762a1bSJed Brown   PetscFunctionReturn(0);
914*c4762a1bSJed Brown }
915*c4762a1bSJed Brown 
916*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C2(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
917*c4762a1bSJed Brown {
918*c4762a1bSJed Brown   PetscErrorCode    ierr;
919*c4762a1bSJed Brown   const PetscScalar *y;
920*c4762a1bSJed Brown   PetscInt          N,i,col[2];
921*c4762a1bSJed Brown   PetscScalar       value[2];
922*c4762a1bSJed Brown 
923*c4762a1bSJed Brown   PetscFunctionBegin;
924*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
925*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
926*c4762a1bSJed Brown   i = 0;
927*c4762a1bSJed Brown   value[0] = a+1;                 col[0] = 0;
928*c4762a1bSJed Brown   value[1] = 0;                   col[1] = 1;
929*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
930*c4762a1bSJed Brown   for (i = 0; i < N; i++) {
931*c4762a1bSJed Brown     value[0] = -(PetscReal) i;      col[0] = i-1;
932*c4762a1bSJed Brown     value[1] = a+(PetscReal)(i+1);  col[1] = i;
933*c4762a1bSJed Brown     ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
934*c4762a1bSJed Brown   }
935*c4762a1bSJed Brown   i = N-1;
936*c4762a1bSJed Brown   value[0] = -(PetscReal) (N-1);  col[0] = N-2;
937*c4762a1bSJed Brown   value[1] = a;                   col[1] = N-1;
938*c4762a1bSJed Brown   ierr = MatSetValues(A,1,&i,2,col,value,INSERT_VALUES);CHKERRQ(ierr);
939*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
940*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
941*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
942*c4762a1bSJed Brown   PetscFunctionReturn(0);
943*c4762a1bSJed Brown }
944*c4762a1bSJed Brown 
945*c4762a1bSJed Brown /* Hull, 1972, Problem C3 and C4 */
946*c4762a1bSJed Brown 
947*c4762a1bSJed Brown PetscErrorCode RHSFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec F, void *s)
948*c4762a1bSJed Brown {
949*c4762a1bSJed Brown   PetscErrorCode    ierr;
950*c4762a1bSJed Brown   PetscScalar       *f;
951*c4762a1bSJed Brown   const PetscScalar *y;
952*c4762a1bSJed Brown   PetscInt          N,i;
953*c4762a1bSJed Brown 
954*c4762a1bSJed Brown   PetscFunctionBegin;
955*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
956*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
957*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
958*c4762a1bSJed Brown   f[0] = -2.0*y[0] + y[1];
959*c4762a1bSJed Brown   for (i = 1; i < N-1; i++) {
960*c4762a1bSJed Brown     f[i] = y[i-1] - 2.0*y[i] + y[i+1];
961*c4762a1bSJed Brown   }
962*c4762a1bSJed Brown   f[N-1] = y[N-2] - 2.0*y[N-1];
963*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
964*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
965*c4762a1bSJed Brown   PetscFunctionReturn(0);
966*c4762a1bSJed Brown }
967*c4762a1bSJed Brown 
968*c4762a1bSJed Brown PetscErrorCode IFunction_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, Vec F, void *s)
969*c4762a1bSJed Brown {
970*c4762a1bSJed Brown   PetscErrorCode    ierr;
971*c4762a1bSJed Brown   PetscScalar       *f;
972*c4762a1bSJed Brown   const PetscScalar *y;
973*c4762a1bSJed Brown   PetscInt          N,i;
974*c4762a1bSJed Brown 
975*c4762a1bSJed Brown   PetscFunctionBegin;
976*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
977*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
978*c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
979*c4762a1bSJed Brown   f[0] = -2.0*y[0] + y[1];
980*c4762a1bSJed Brown   for (i = 1; i < N-1; i++) {
981*c4762a1bSJed Brown     f[i] = y[i-1] - 2.0*y[i] + y[i+1];
982*c4762a1bSJed Brown   }
983*c4762a1bSJed Brown   f[N-1] = y[N-2] - 2.0*y[N-1];
984*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
985*c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
986*c4762a1bSJed Brown   /* Left hand side = ydot - f(y) */
987*c4762a1bSJed Brown   ierr = VecAYPX(F,-1.0,Ydot);CHKERRQ(ierr);
988*c4762a1bSJed Brown   PetscFunctionReturn(0);
989*c4762a1bSJed Brown }
990*c4762a1bSJed Brown 
991*c4762a1bSJed Brown PetscErrorCode IJacobian_Hull1972C34(TS ts, PetscReal t, Vec Y, Vec Ydot, PetscReal a, Mat A, Mat B, void *s)
992*c4762a1bSJed Brown {
993*c4762a1bSJed Brown   PetscErrorCode    ierr;
994*c4762a1bSJed Brown   const PetscScalar *y;
995*c4762a1bSJed Brown   PetscScalar       value[3];
996*c4762a1bSJed Brown   PetscInt          N,i,col[3];
997*c4762a1bSJed Brown 
998*c4762a1bSJed Brown   PetscFunctionBegin;
999*c4762a1bSJed Brown   ierr = VecGetSize (Y,&N);CHKERRQ(ierr);
1000*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y,&y);CHKERRQ(ierr);
1001*c4762a1bSJed Brown   for (i = 0; i < N; i++) {
1002*c4762a1bSJed Brown     if (i == 0) {
1003*c4762a1bSJed Brown       value[0] = a+2;  col[0] = i;
1004*c4762a1bSJed Brown       value[1] =  -1;  col[1] = i+1;
1005*c4762a1bSJed Brown       value[2] =  0;   col[2] = i+2;
1006*c4762a1bSJed Brown     } else if (i == N-1) {
1007*c4762a1bSJed Brown       value[0] =  0;   col[0] = i-2;
1008*c4762a1bSJed Brown       value[1] =  -1;  col[1] = i-1;
1009*c4762a1bSJed Brown       value[2] = a+2;  col[2] = i;
1010*c4762a1bSJed Brown     } else {
1011*c4762a1bSJed Brown       value[0] = -1;   col[0] = i-1;
1012*c4762a1bSJed Brown       value[1] = a+2;  col[1] = i;
1013*c4762a1bSJed Brown       value[2] = -1;   col[2] = i+1;
1014*c4762a1bSJed Brown     }
1015*c4762a1bSJed Brown     ierr = MatSetValues(A,1,&i,3,col,value,INSERT_VALUES);CHKERRQ(ierr);
1016*c4762a1bSJed Brown   }
1017*c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1018*c4762a1bSJed Brown   ierr = MatAssemblyEnd  (A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1019*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y,&y);CHKERRQ(ierr);
1020*c4762a1bSJed Brown   PetscFunctionReturn(0);
1021*c4762a1bSJed Brown }
1022*c4762a1bSJed Brown 
1023*c4762a1bSJed Brown /***************************************************************************/
1024*c4762a1bSJed Brown 
1025*c4762a1bSJed Brown /* Sets the initial solution for the IVP and sets up the function pointers*/
1026*c4762a1bSJed Brown PetscErrorCode Initialize(Vec Y, void* s)
1027*c4762a1bSJed Brown {
1028*c4762a1bSJed Brown   PetscErrorCode ierr;
1029*c4762a1bSJed Brown   char          *p = (char*) s;
1030*c4762a1bSJed Brown   PetscScalar   *y;
1031*c4762a1bSJed Brown   PetscReal     t0;
1032*c4762a1bSJed Brown   PetscInt      N = GetSize((const char *)s);
1033*c4762a1bSJed Brown   PetscBool     flg;
1034*c4762a1bSJed Brown 
1035*c4762a1bSJed Brown   PetscFunctionBegin;
1036*c4762a1bSJed Brown   VecZeroEntries(Y);
1037*c4762a1bSJed Brown   ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1038*c4762a1bSJed Brown   if (!strcmp(p,"hull1972a1")) {
1039*c4762a1bSJed Brown     y[0] = 1.0;
1040*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A1;
1041*c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A1;
1042*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A1;
1043*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A1;
1044*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972a2")) {
1045*c4762a1bSJed Brown     y[0] = 1.0;
1046*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A2;
1047*c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A2;
1048*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A2;
1049*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A2;
1050*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972a3")) {
1051*c4762a1bSJed Brown     y[0] = 1.0;
1052*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A3;
1053*c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A3;
1054*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A3;
1055*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A3;
1056*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972a4")) {
1057*c4762a1bSJed Brown     y[0] = 1.0;
1058*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A4;
1059*c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A4;
1060*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A4;
1061*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A4;
1062*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972a5")) {
1063*c4762a1bSJed Brown     y[0] = 4.0;
1064*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972A5;
1065*c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Hull1972A5;
1066*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972A5;
1067*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972A5;
1068*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972b1")) {
1069*c4762a1bSJed Brown     y[0] = 1.0;
1070*c4762a1bSJed Brown     y[1] = 3.0;
1071*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B1;
1072*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B1;
1073*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B1;
1074*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972b2")) {
1075*c4762a1bSJed Brown     y[0] = 2.0;
1076*c4762a1bSJed Brown     y[1] = 0.0;
1077*c4762a1bSJed Brown     y[2] = 1.0;
1078*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B2;
1079*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B2;
1080*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B2;
1081*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972b3")) {
1082*c4762a1bSJed Brown     y[0] = 1.0;
1083*c4762a1bSJed Brown     y[1] = 0.0;
1084*c4762a1bSJed Brown     y[2] = 0.0;
1085*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B3;
1086*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B3;
1087*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B3;
1088*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972b4")) {
1089*c4762a1bSJed Brown     y[0] = 3.0;
1090*c4762a1bSJed Brown     y[1] = 0.0;
1091*c4762a1bSJed Brown     y[2] = 0.0;
1092*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B4;
1093*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B4;
1094*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B4;
1095*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972b5")) {
1096*c4762a1bSJed Brown     y[0] = 0.0;
1097*c4762a1bSJed Brown     y[1] = 1.0;
1098*c4762a1bSJed Brown     y[2] = 1.0;
1099*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972B5;
1100*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972B5;
1101*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972B5;
1102*c4762a1bSJed Brown   } else if (!strcmp(p,"kulik2013i")) {
1103*c4762a1bSJed Brown     t0=0.;
1104*c4762a1bSJed Brown     y[0] = PetscExpReal(PetscSinReal(t0*t0));
1105*c4762a1bSJed Brown     y[1] = PetscExpReal(5.*PetscSinReal(t0*t0));
1106*c4762a1bSJed Brown     y[2] = PetscSinReal(t0*t0)+1.0;
1107*c4762a1bSJed Brown     y[3] = PetscCosReal(t0*t0);
1108*c4762a1bSJed Brown     RHSFunction = RHSFunction_Kulikov2013I;
1109*c4762a1bSJed Brown     RHSJacobian = RHSJacobian_Kulikov2013I;
1110*c4762a1bSJed Brown     IFunction   = IFunction_Kulikov2013I;
1111*c4762a1bSJed Brown     IJacobian   = IJacobian_Kulikov2013I;
1112*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972c1")) {
1113*c4762a1bSJed Brown     y[0] = 1.0;
1114*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C1;
1115*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C1;
1116*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C1;
1117*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972c2")) {
1118*c4762a1bSJed Brown     y[0] = 1.0;
1119*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C2;
1120*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C2;
1121*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C2;
1122*c4762a1bSJed Brown   } else if ((!strcmp(p,"hull1972c3"))
1123*c4762a1bSJed Brown            ||(!strcmp(p,"hull1972c4"))){
1124*c4762a1bSJed Brown     y[0] = 1.0;
1125*c4762a1bSJed Brown     RHSFunction = RHSFunction_Hull1972C34;
1126*c4762a1bSJed Brown     IFunction   = IFunction_Hull1972C34;
1127*c4762a1bSJed Brown     IJacobian   = IJacobian_Hull1972C34;
1128*c4762a1bSJed Brown   }
1129*c4762a1bSJed Brown   ierr = PetscOptionsGetScalarArray(NULL,NULL,"-yinit",y,&N,&flg);CHKERRQ(ierr);
1130*c4762a1bSJed Brown   if ((N != GetSize((const char*)s)) && flg) SETERRQ2(PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Number of initial values %D does not match problem size %D.\n",N,GetSize((const char*)s));
1131*c4762a1bSJed Brown   ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1132*c4762a1bSJed Brown   PetscFunctionReturn(0);
1133*c4762a1bSJed Brown }
1134*c4762a1bSJed Brown 
1135*c4762a1bSJed Brown /* Calculates the exact solution to problems that have one */
1136*c4762a1bSJed Brown PetscErrorCode ExactSolution(Vec Y, void* s, PetscReal t, PetscBool *flag)
1137*c4762a1bSJed Brown {
1138*c4762a1bSJed Brown   PetscErrorCode ierr;
1139*c4762a1bSJed Brown   char          *p = (char*) s;
1140*c4762a1bSJed Brown   PetscScalar   *y;
1141*c4762a1bSJed Brown 
1142*c4762a1bSJed Brown   PetscFunctionBegin;
1143*c4762a1bSJed Brown   if (!strcmp(p,"hull1972a1")) {
1144*c4762a1bSJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1145*c4762a1bSJed Brown     y[0] = PetscExpReal(-t);
1146*c4762a1bSJed Brown     *flag = PETSC_TRUE;
1147*c4762a1bSJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1148*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972a2")) {
1149*c4762a1bSJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1150*c4762a1bSJed Brown     y[0] = 1.0/PetscSqrtReal(t+1);
1151*c4762a1bSJed Brown     *flag = PETSC_TRUE;
1152*c4762a1bSJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1153*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972a3")) {
1154*c4762a1bSJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1155*c4762a1bSJed Brown     y[0] = PetscExpReal(PetscSinReal(t));
1156*c4762a1bSJed Brown     *flag = PETSC_TRUE;
1157*c4762a1bSJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1158*c4762a1bSJed Brown   } else if (!strcmp(p,"hull1972a4")) {
1159*c4762a1bSJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1160*c4762a1bSJed Brown     y[0] = 20.0/(1+19.0*PetscExpReal(-t/4.0));
1161*c4762a1bSJed Brown     *flag = PETSC_TRUE;
1162*c4762a1bSJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1163*c4762a1bSJed Brown   } else if (!strcmp(p,"kulik2013i")) {
1164*c4762a1bSJed Brown     ierr = VecGetArray(Y,&y);CHKERRQ(ierr);
1165*c4762a1bSJed Brown     y[0] = PetscExpReal(PetscSinReal(t*t));
1166*c4762a1bSJed Brown     y[1] = PetscExpReal(5.*PetscSinReal(t*t));
1167*c4762a1bSJed Brown     y[2] = PetscSinReal(t*t)+1.0;
1168*c4762a1bSJed Brown     y[3] = PetscCosReal(t*t);
1169*c4762a1bSJed Brown     *flag = PETSC_TRUE;
1170*c4762a1bSJed Brown     ierr = VecRestoreArray(Y,&y);CHKERRQ(ierr);
1171*c4762a1bSJed Brown   } else {
1172*c4762a1bSJed Brown     ierr = VecSet(Y,0);CHKERRQ(ierr);
1173*c4762a1bSJed Brown     *flag = PETSC_FALSE;
1174*c4762a1bSJed Brown   }
1175*c4762a1bSJed Brown   PetscFunctionReturn(0);
1176*c4762a1bSJed Brown }
1177*c4762a1bSJed Brown 
1178*c4762a1bSJed Brown /* Solves the specified ODE and computes the error if exact solution is available */
1179*c4762a1bSJed Brown PetscErrorCode SolveODE(char* ptype, PetscReal dt, PetscReal tfinal, PetscInt maxiter, PetscReal *error, PetscBool *exact_flag)
1180*c4762a1bSJed Brown {
1181*c4762a1bSJed Brown   PetscErrorCode  ierr;             /* Error code                             */
1182*c4762a1bSJed Brown   TS              ts;               /* time-integrator                        */
1183*c4762a1bSJed Brown   Vec             Y;                /* Solution vector                        */
1184*c4762a1bSJed Brown   Vec             Yex;              /* Exact solution                         */
1185*c4762a1bSJed Brown   PetscInt        N;                /* Size of the system of equations        */
1186*c4762a1bSJed Brown   TSType          time_scheme;      /* Type of time-integration scheme        */
1187*c4762a1bSJed Brown   Mat             Jac = NULL;       /* Jacobian matrix                        */
1188*c4762a1bSJed Brown   Vec             Yerr;             /* Auxiliary solution vector              */
1189*c4762a1bSJed Brown   PetscReal       err_norm;         /* Estimated error norm                   */
1190*c4762a1bSJed Brown   PetscReal       final_time;       /* Actual final time from the integrator  */
1191*c4762a1bSJed Brown 
1192*c4762a1bSJed Brown   PetscFunctionBegin;
1193*c4762a1bSJed Brown   N = GetSize((const char *)&ptype[0]);
1194*c4762a1bSJed Brown   if (N < 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_ARG_SIZ,"Illegal problem specification.\n");
1195*c4762a1bSJed Brown   ierr = VecCreate(PETSC_COMM_WORLD,&Y);CHKERRQ(ierr);
1196*c4762a1bSJed Brown   ierr = VecSetSizes(Y,N,PETSC_DECIDE);CHKERRQ(ierr);
1197*c4762a1bSJed Brown   ierr = VecSetUp(Y);CHKERRQ(ierr);
1198*c4762a1bSJed Brown   ierr = VecSet(Y,0);CHKERRQ(ierr);
1199*c4762a1bSJed Brown 
1200*c4762a1bSJed Brown   /* Initialize the problem */
1201*c4762a1bSJed Brown   ierr = Initialize(Y,&ptype[0]);CHKERRQ(ierr);
1202*c4762a1bSJed Brown 
1203*c4762a1bSJed Brown   /* Create and initialize the time-integrator                            */
1204*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
1205*c4762a1bSJed Brown   /* Default time integration options                                     */
1206*c4762a1bSJed Brown   ierr = TSSetType(ts,TSRK);CHKERRQ(ierr);
1207*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts,maxiter);CHKERRQ(ierr);
1208*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,tfinal);CHKERRQ(ierr);
1209*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,dt);CHKERRQ(ierr);
1210*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
1211*c4762a1bSJed Brown   /* Read command line options for time integration                       */
1212*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
1213*c4762a1bSJed Brown   /* Set solution vector                                                  */
1214*c4762a1bSJed Brown   ierr = TSSetSolution(ts,Y);CHKERRQ(ierr);
1215*c4762a1bSJed Brown   /* Specify left/right-hand side functions                               */
1216*c4762a1bSJed Brown   ierr = TSGetType(ts,&time_scheme);CHKERRQ(ierr);
1217*c4762a1bSJed Brown 
1218*c4762a1bSJed Brown   if ((!strcmp(time_scheme,TSEULER)) || (!strcmp(time_scheme,TSRK)) || (!strcmp(time_scheme,TSSSP) || (!strcmp(time_scheme,TSGLEE)))) {
1219*c4762a1bSJed Brown     /* Explicit time-integration -> specify right-hand side function ydot = f(y) */
1220*c4762a1bSJed Brown     ierr = TSSetRHSFunction(ts,NULL,RHSFunction,&ptype[0]);CHKERRQ(ierr);
1221*c4762a1bSJed Brown     ierr = MatCreate(PETSC_COMM_WORLD,&Jac);CHKERRQ(ierr);
1222*c4762a1bSJed Brown     ierr = MatSetSizes(Jac,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
1223*c4762a1bSJed Brown     ierr = MatSetFromOptions(Jac);CHKERRQ(ierr);
1224*c4762a1bSJed Brown     ierr = MatSetUp(Jac);CHKERRQ(ierr);
1225*c4762a1bSJed Brown     ierr = TSSetRHSJacobian(ts,Jac,Jac,RHSJacobian,&ptype[0]);CHKERRQ(ierr);
1226*c4762a1bSJed Brown   } else if ((!strcmp(time_scheme,TSTHETA)) || (!strcmp(time_scheme,TSBEULER)) || (!strcmp(time_scheme,TSCN)) || (!strcmp(time_scheme,TSALPHA)) || (!strcmp(time_scheme,TSARKIMEX))) {
1227*c4762a1bSJed Brown     /* Implicit time-integration -> specify left-hand side function ydot-f(y) = 0 */
1228*c4762a1bSJed Brown     /* and its Jacobian function                                                 */
1229*c4762a1bSJed Brown     ierr = TSSetIFunction(ts,NULL,IFunction,&ptype[0]);CHKERRQ(ierr);
1230*c4762a1bSJed Brown     ierr = MatCreate(PETSC_COMM_WORLD,&Jac);CHKERRQ(ierr);
1231*c4762a1bSJed Brown     ierr = MatSetSizes(Jac,PETSC_DECIDE,PETSC_DECIDE,N,N);CHKERRQ(ierr);
1232*c4762a1bSJed Brown     ierr = MatSetFromOptions(Jac);CHKERRQ(ierr);
1233*c4762a1bSJed Brown     ierr = MatSetUp(Jac);CHKERRQ(ierr);
1234*c4762a1bSJed Brown     ierr = TSSetIJacobian(ts,Jac,Jac,IJacobian,&ptype[0]);CHKERRQ(ierr);
1235*c4762a1bSJed Brown   }
1236*c4762a1bSJed Brown 
1237*c4762a1bSJed Brown   /* Solve */
1238*c4762a1bSJed Brown   ierr = TSSolve(ts,Y);CHKERRQ(ierr);
1239*c4762a1bSJed Brown   ierr = TSGetTime(ts,&final_time);CHKERRQ(ierr);
1240*c4762a1bSJed Brown 
1241*c4762a1bSJed Brown   /* Get the estimated error, if available */
1242*c4762a1bSJed Brown   ierr = VecDuplicate(Y,&Yerr);CHKERRQ(ierr);
1243*c4762a1bSJed Brown   ierr = VecZeroEntries(Yerr);CHKERRQ(ierr);
1244*c4762a1bSJed Brown   ierr = TSGetTimeError(ts,0,&Yerr);CHKERRQ(ierr);
1245*c4762a1bSJed Brown   ierr = VecNorm(Yerr,NORM_2,&err_norm);CHKERRQ(ierr);
1246*c4762a1bSJed Brown   ierr = VecDestroy(&Yerr);CHKERRQ(ierr);
1247*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD,"Estimated Error = %E.\n",err_norm);CHKERRQ(ierr);
1248*c4762a1bSJed Brown 
1249*c4762a1bSJed Brown   /* Exact solution */
1250*c4762a1bSJed Brown   ierr = VecDuplicate(Y,&Yex);CHKERRQ(ierr);
1251*c4762a1bSJed Brown   if(PetscAbsScalar(final_time-tfinal)>2.*PETSC_MACHINE_EPSILON) {
1252*c4762a1bSJed Brown     ierr = 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);CHKERRQ(ierr);
1253*c4762a1bSJed Brown   }
1254*c4762a1bSJed Brown   ierr = ExactSolution(Yex,&ptype[0],final_time,exact_flag);CHKERRQ(ierr);
1255*c4762a1bSJed Brown 
1256*c4762a1bSJed Brown   /* Calculate Error */
1257*c4762a1bSJed Brown   ierr = VecAYPX(Yex,-1.0,Y);CHKERRQ(ierr);
1258*c4762a1bSJed Brown   ierr = VecNorm(Yex,NORM_2,error);CHKERRQ(ierr);
1259*c4762a1bSJed Brown   *error = PetscSqrtReal(((*error)*(*error))/N);
1260*c4762a1bSJed Brown 
1261*c4762a1bSJed Brown   /* Clean up and finalize */
1262*c4762a1bSJed Brown   ierr = MatDestroy(&Jac);CHKERRQ(ierr);
1263*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
1264*c4762a1bSJed Brown   ierr = VecDestroy(&Yex);CHKERRQ(ierr);
1265*c4762a1bSJed Brown   ierr = VecDestroy(&Y);CHKERRQ(ierr);
1266*c4762a1bSJed Brown 
1267*c4762a1bSJed Brown   PetscFunctionReturn(0);
1268*c4762a1bSJed Brown }
1269*c4762a1bSJed Brown 
1270*c4762a1bSJed Brown int main(int argc, char **argv)
1271*c4762a1bSJed Brown {
1272*c4762a1bSJed Brown   PetscErrorCode  ierr;                       /* Error code                                           */
1273*c4762a1bSJed Brown   char            ptype[256] = "hull1972a1";  /* Problem specification                                */
1274*c4762a1bSJed Brown   PetscInt        n_refine   = 1;             /* Number of refinement levels for convergence analysis */
1275*c4762a1bSJed Brown   PetscReal       refine_fac = 2.0;           /* Refinement factor for dt                             */
1276*c4762a1bSJed Brown   PetscReal       dt_initial = 0.01;          /* Initial default value of dt                          */
1277*c4762a1bSJed Brown   PetscReal       dt;
1278*c4762a1bSJed Brown   PetscReal       tfinal     = 20.0;          /* Final time for the time-integration                  */
1279*c4762a1bSJed Brown   PetscInt        maxiter    = 100000;        /* Maximum number of time-integration iterations        */
1280*c4762a1bSJed Brown   PetscReal       *error;                     /* Array to store the errors for convergence analysis   */
1281*c4762a1bSJed Brown   PetscMPIInt     size;                      /* No of processors                                     */
1282*c4762a1bSJed Brown   PetscBool       flag;                       /* Flag denoting availability of exact solution         */
1283*c4762a1bSJed Brown   PetscInt        r;
1284*c4762a1bSJed Brown 
1285*c4762a1bSJed Brown   /* Initialize program */
1286*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
1287*c4762a1bSJed Brown 
1288*c4762a1bSJed Brown   /* Check if running with only 1 proc */
1289*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
1290*c4762a1bSJed Brown   if (size>1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");
1291*c4762a1bSJed Brown 
1292*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,NULL,NULL);CHKERRQ(ierr);
1293*c4762a1bSJed Brown   ierr = PetscOptionsString("-problem","Problem specification","<hull1972a1>",ptype,ptype,sizeof(ptype),NULL);CHKERRQ(ierr);
1294*c4762a1bSJed Brown   ierr = PetscOptionsInt("-refinement_levels","Number of refinement levels for convergence analysis","<1>",n_refine,&n_refine,NULL);CHKERRQ(ierr);
1295*c4762a1bSJed Brown   ierr = PetscOptionsReal("-refinement_factor","Refinement factor for dt","<2.0>",refine_fac,&refine_fac,NULL);CHKERRQ(ierr);
1296*c4762a1bSJed Brown   ierr = PetscOptionsReal("-dt","Time step size (for convergence analysis, initial time step)","<0.01>",dt_initial,&dt_initial,NULL);CHKERRQ(ierr);
1297*c4762a1bSJed Brown   ierr = PetscOptionsReal("-final_time","Final time for the time-integration","<20.0>",tfinal,&tfinal,NULL);CHKERRQ(ierr);
1298*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
1299*c4762a1bSJed Brown 
1300*c4762a1bSJed Brown   ierr = PetscMalloc1(n_refine,&error);CHKERRQ(ierr);
1301*c4762a1bSJed Brown   for (r = 0,dt = dt_initial; r < n_refine; r++) {
1302*c4762a1bSJed Brown     error[r] = 0;
1303*c4762a1bSJed Brown     if (r > 0) dt /= refine_fac;
1304*c4762a1bSJed Brown 
1305*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"Solving ODE \"%s\" with dt %f, final time %f and system size %D.\n",ptype,(double)dt,(double)tfinal,GetSize(&ptype[0]));CHKERRQ(ierr);
1306*c4762a1bSJed Brown     ierr = SolveODE(&ptype[0],dt,tfinal,maxiter,&error[r],&flag);CHKERRQ(ierr);
1307*c4762a1bSJed Brown     if (flag) {
1308*c4762a1bSJed Brown       /* If exact solution available for the specified ODE */
1309*c4762a1bSJed Brown       if (r > 0) {
1310*c4762a1bSJed Brown         PetscReal conv_rate = (PetscLogReal(error[r]) - PetscLogReal(error[r-1])) / (-PetscLogReal(refine_fac));
1311*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_WORLD,"Error           = %E,\tConvergence rate = %f.\n",(double)error[r],(double)conv_rate);CHKERRQ(ierr);
1312*c4762a1bSJed Brown       } else {
1313*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_WORLD,"Error           = %E.\n",error[r]);CHKERRQ(ierr);
1314*c4762a1bSJed Brown       }
1315*c4762a1bSJed Brown     }
1316*c4762a1bSJed Brown   }
1317*c4762a1bSJed Brown   ierr = PetscFree(error);CHKERRQ(ierr);
1318*c4762a1bSJed Brown   ierr = PetscFinalize();
1319*c4762a1bSJed Brown   return ierr;
1320*c4762a1bSJed Brown }
1321*c4762a1bSJed Brown 
1322*c4762a1bSJed Brown /*TEST
1323*c4762a1bSJed Brown 
1324*c4762a1bSJed Brown     test:
1325*c4762a1bSJed Brown       suffix: 2
1326*c4762a1bSJed Brown       args: -ts_type glee -final_time 5 -ts_adapt_type none
1327*c4762a1bSJed Brown       timeoutfactor: 3
1328*c4762a1bSJed Brown       requires: !single
1329*c4762a1bSJed Brown 
1330*c4762a1bSJed Brown     test:
1331*c4762a1bSJed Brown       suffix: 3
1332*c4762a1bSJed 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
1333*c4762a1bSJed Brown       timeoutfactor: 3
1334*c4762a1bSJed Brown       requires: !single
1335*c4762a1bSJed Brown 
1336*c4762a1bSJed Brown     test:
1337*c4762a1bSJed Brown       suffix: 4
1338*c4762a1bSJed 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
1339*c4762a1bSJed Brown       timeoutfactor: 3
1340*c4762a1bSJed Brown       requires: !single !__float128
1341*c4762a1bSJed Brown 
1342*c4762a1bSJed Brown TEST*/
1343