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