xref: /petsc/src/ts/tutorials/ex8.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Nonlinear DAE benchmark problems.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*
5c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
6c4762a1bSJed Brown    file automatically includes:
7c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
8c4762a1bSJed Brown      petscmat.h - matrices
9c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
10c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
11c4762a1bSJed Brown      petscksp.h   - linear solvers
12c4762a1bSJed Brown */
13c4762a1bSJed Brown #include <petscts.h>
14c4762a1bSJed Brown 
15c4762a1bSJed Brown typedef struct _Problem *Problem;
16c4762a1bSJed Brown struct _Problem {
17c4762a1bSJed Brown   PetscErrorCode (*destroy)(Problem);
18c4762a1bSJed Brown   TSIFunction function;
19c4762a1bSJed Brown   TSIJacobian jacobian;
20c4762a1bSJed Brown   PetscErrorCode (*solution)(PetscReal, Vec, void *);
21c4762a1bSJed Brown   MPI_Comm  comm;
22c4762a1bSJed Brown   PetscReal final_time;
23c4762a1bSJed Brown   PetscInt  n;
24c4762a1bSJed Brown   PetscBool hasexact;
25c4762a1bSJed Brown   void     *data;
26c4762a1bSJed Brown };
27c4762a1bSJed Brown 
28c4762a1bSJed Brown /*
29c4762a1bSJed Brown       Stiff 3-variable system from chemical reactions, due to Robertson (1966), problem ROBER in Hairer&Wanner, ODE 2, 1996
30c4762a1bSJed Brown */
319371c9d4SSatish Balay static PetscErrorCode RoberFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
32c4762a1bSJed Brown   PetscScalar       *f;
33c4762a1bSJed Brown   const PetscScalar *x, *xdot;
34c4762a1bSJed Brown 
35c4762a1bSJed Brown   PetscFunctionBeginUser;
369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
389566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
39c4762a1bSJed Brown   f[0] = xdot[0] + 0.04 * x[0] - 1e4 * x[1] * x[2];
40c4762a1bSJed Brown   f[1] = xdot[1] - 0.04 * x[0] + 1e4 * x[1] * x[2] + 3e7 * PetscSqr(x[1]);
41c4762a1bSJed Brown   f[2] = xdot[2] - 3e7 * PetscSqr(x[1]);
429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
45c4762a1bSJed Brown   PetscFunctionReturn(0);
46c4762a1bSJed Brown }
47c4762a1bSJed Brown 
489371c9d4SSatish Balay static PetscErrorCode RoberJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
49c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1, 2};
50c4762a1bSJed Brown   PetscScalar        J[3][3];
51c4762a1bSJed Brown   const PetscScalar *x, *xdot;
52c4762a1bSJed Brown 
53c4762a1bSJed Brown   PetscFunctionBeginUser;
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
569371c9d4SSatish Balay   J[0][0] = a + 0.04;
579371c9d4SSatish Balay   J[0][1] = -1e4 * x[2];
589371c9d4SSatish Balay   J[0][2] = -1e4 * x[1];
599371c9d4SSatish Balay   J[1][0] = -0.04;
609371c9d4SSatish Balay   J[1][1] = a + 1e4 * x[2] + 3e7 * 2 * x[1];
619371c9d4SSatish Balay   J[1][2] = 1e4 * x[1];
629371c9d4SSatish Balay   J[2][0] = 0;
639371c9d4SSatish Balay   J[2][1] = -3e7 * 2 * x[1];
649371c9d4SSatish Balay   J[2][2] = a;
659566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 3, rowcol, 3, rowcol, &J[0][0], INSERT_VALUES));
669566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
68c4762a1bSJed Brown 
699566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
709566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
71c4762a1bSJed Brown   if (A != B) {
729566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
739566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
74c4762a1bSJed Brown   }
75c4762a1bSJed Brown   PetscFunctionReturn(0);
76c4762a1bSJed Brown }
77c4762a1bSJed Brown 
789371c9d4SSatish Balay static PetscErrorCode RoberSolution(PetscReal t, Vec X, void *ctx) {
79c4762a1bSJed Brown   PetscScalar *x;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   PetscFunctionBeginUser;
823c633725SBarry Smith   PetscCheck(t == 0, PETSC_COMM_WORLD, PETSC_ERR_SUP, "not implemented");
839566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
84c4762a1bSJed Brown   x[0] = 1;
85c4762a1bSJed Brown   x[1] = 0;
86c4762a1bSJed Brown   x[2] = 0;
879566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
88c4762a1bSJed Brown   PetscFunctionReturn(0);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
919371c9d4SSatish Balay static PetscErrorCode RoberCreate(Problem p) {
92c4762a1bSJed Brown   PetscFunctionBeginUser;
93c4762a1bSJed Brown   p->destroy    = 0;
94c4762a1bSJed Brown   p->function   = &RoberFunction;
95c4762a1bSJed Brown   p->jacobian   = &RoberJacobian;
96c4762a1bSJed Brown   p->solution   = &RoberSolution;
97c4762a1bSJed Brown   p->final_time = 1e11;
98c4762a1bSJed Brown   p->n          = 3;
99c4762a1bSJed Brown   PetscFunctionReturn(0);
100c4762a1bSJed Brown }
101c4762a1bSJed Brown 
102c4762a1bSJed Brown /*
103c4762a1bSJed Brown      Stiff scalar valued problem
104c4762a1bSJed Brown */
105c4762a1bSJed Brown 
106c4762a1bSJed Brown typedef struct {
107c4762a1bSJed Brown   PetscReal lambda;
108c4762a1bSJed Brown } CECtx;
109c4762a1bSJed Brown 
1109371c9d4SSatish Balay static PetscErrorCode CEDestroy(Problem p) {
111c4762a1bSJed Brown   PetscFunctionBeginUser;
1129566063dSJacob Faibussowitsch   PetscCall(PetscFree(p->data));
113c4762a1bSJed Brown   PetscFunctionReturn(0);
114c4762a1bSJed Brown }
115c4762a1bSJed Brown 
1169371c9d4SSatish Balay static PetscErrorCode CEFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
117c4762a1bSJed Brown   PetscReal          l = ((CECtx *)ctx)->lambda;
118c4762a1bSJed Brown   PetscScalar       *f;
119c4762a1bSJed Brown   const PetscScalar *x, *xdot;
120c4762a1bSJed Brown 
121c4762a1bSJed Brown   PetscFunctionBeginUser;
1229566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1239566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
1249566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
125c4762a1bSJed Brown   f[0] = xdot[0] + l * (x[0] - PetscCosReal(t));
126c4762a1bSJed Brown #if 0
1279566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD," f(t=%g,x=%g,xdot=%g) = %g\n",(double)t,(double)x[0],(double)xdot[0],(double)f[0]));
128c4762a1bSJed Brown #endif
1299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
1319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
132c4762a1bSJed Brown   PetscFunctionReturn(0);
133c4762a1bSJed Brown }
134c4762a1bSJed Brown 
1359371c9d4SSatish Balay static PetscErrorCode CEJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
136c4762a1bSJed Brown   PetscReal          l        = ((CECtx *)ctx)->lambda;
137c4762a1bSJed Brown   PetscInt           rowcol[] = {0};
138c4762a1bSJed Brown   PetscScalar        J[1][1];
139c4762a1bSJed Brown   const PetscScalar *x, *xdot;
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   PetscFunctionBeginUser;
1429566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
144c4762a1bSJed Brown   J[0][0] = a + l;
1459566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, rowcol, 1, rowcol, &J[0][0], INSERT_VALUES));
1469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
148c4762a1bSJed Brown 
1499566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
151c4762a1bSJed Brown   if (A != B) {
1529566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
154c4762a1bSJed Brown   }
155c4762a1bSJed Brown   PetscFunctionReturn(0);
156c4762a1bSJed Brown }
157c4762a1bSJed Brown 
1589371c9d4SSatish Balay static PetscErrorCode CESolution(PetscReal t, Vec X, void *ctx) {
159c4762a1bSJed Brown   PetscReal    l = ((CECtx *)ctx)->lambda;
160c4762a1bSJed Brown   PetscScalar *x;
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   PetscFunctionBeginUser;
1639566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
164c4762a1bSJed Brown   x[0] = l / (l * l + 1) * (l * PetscCosReal(t) + PetscSinReal(t)) - l * l / (l * l + 1) * PetscExpReal(-l * t);
1659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
166c4762a1bSJed Brown   PetscFunctionReturn(0);
167c4762a1bSJed Brown }
168c4762a1bSJed Brown 
1699371c9d4SSatish Balay static PetscErrorCode CECreate(Problem p) {
170c4762a1bSJed Brown   CECtx *ce;
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   PetscFunctionBeginUser;
1739566063dSJacob Faibussowitsch   PetscCall(PetscMalloc(sizeof(CECtx), &ce));
174c4762a1bSJed Brown   p->data = (void *)ce;
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   p->destroy    = &CEDestroy;
177c4762a1bSJed Brown   p->function   = &CEFunction;
178c4762a1bSJed Brown   p->jacobian   = &CEJacobian;
179c4762a1bSJed Brown   p->solution   = &CESolution;
180c4762a1bSJed Brown   p->final_time = 10;
181c4762a1bSJed Brown   p->n          = 1;
182c4762a1bSJed Brown   p->hasexact   = PETSC_TRUE;
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   ce->lambda = 10;
185d0609cedSBarry Smith   PetscOptionsBegin(p->comm, NULL, "CE options", "");
1869371c9d4SSatish Balay   { PetscCall(PetscOptionsReal("-problem_ce_lambda", "Parameter controlling stiffness: xdot + lambda*(x - cos(t))", "", ce->lambda, &ce->lambda, NULL)); }
187d0609cedSBarry Smith   PetscOptionsEnd();
188c4762a1bSJed Brown   PetscFunctionReturn(0);
189c4762a1bSJed Brown }
190c4762a1bSJed Brown 
191c4762a1bSJed Brown /*
1920e3d61c9SBarry Smith    Stiff 3-variable oscillatory system from chemical reactions. problem OREGO in Hairer&Wanner
193c4762a1bSJed Brown */
1949371c9d4SSatish Balay static PetscErrorCode OregoFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
195c4762a1bSJed Brown   PetscScalar       *f;
196c4762a1bSJed Brown   const PetscScalar *x, *xdot;
197c4762a1bSJed Brown 
198c4762a1bSJed Brown   PetscFunctionBeginUser;
1999566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2009566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
2019566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
202c4762a1bSJed Brown   f[0] = xdot[0] - 77.27 * (x[1] + x[0] * (1. - 8.375e-6 * x[0] - x[1]));
203c4762a1bSJed Brown   f[1] = xdot[1] - 1 / 77.27 * (x[2] - (1. + x[0]) * x[1]);
204c4762a1bSJed Brown   f[2] = xdot[2] - 0.161 * (x[0] - x[2]);
2059566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
2079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
208c4762a1bSJed Brown   PetscFunctionReturn(0);
209c4762a1bSJed Brown }
210c4762a1bSJed Brown 
2119371c9d4SSatish Balay static PetscErrorCode OregoJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
212c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1, 2};
213c4762a1bSJed Brown   PetscScalar        J[3][3];
214c4762a1bSJed Brown   const PetscScalar *x, *xdot;
215c4762a1bSJed Brown 
216c4762a1bSJed Brown   PetscFunctionBeginUser;
2179566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2189566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
219c4762a1bSJed Brown   J[0][0] = a - 77.27 * ((1. - 8.375e-6 * x[0] - x[1]) - 8.375e-6 * x[0]);
220c4762a1bSJed Brown   J[0][1] = -77.27 * (1. - x[0]);
221c4762a1bSJed Brown   J[0][2] = 0;
222c4762a1bSJed Brown   J[1][0] = 1. / 77.27 * x[1];
223c4762a1bSJed Brown   J[1][1] = a + 1. / 77.27 * (1. + x[0]);
224c4762a1bSJed Brown   J[1][2] = -1. / 77.27;
225c4762a1bSJed Brown   J[2][0] = -0.161;
226c4762a1bSJed Brown   J[2][1] = 0;
227c4762a1bSJed Brown   J[2][2] = a + 0.161;
2289566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 3, rowcol, 3, rowcol, &J[0][0], INSERT_VALUES));
2299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
231c4762a1bSJed Brown 
2329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2339566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
234c4762a1bSJed Brown   if (A != B) {
2359566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2369566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
237c4762a1bSJed Brown   }
238c4762a1bSJed Brown   PetscFunctionReturn(0);
239c4762a1bSJed Brown }
240c4762a1bSJed Brown 
2419371c9d4SSatish Balay static PetscErrorCode OregoSolution(PetscReal t, Vec X, void *ctx) {
242c4762a1bSJed Brown   PetscScalar *x;
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   PetscFunctionBeginUser;
2453c633725SBarry Smith   PetscCheck(t == 0, PETSC_COMM_WORLD, PETSC_ERR_SUP, "not implemented");
2469566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
247c4762a1bSJed Brown   x[0] = 1;
248c4762a1bSJed Brown   x[1] = 2;
249c4762a1bSJed Brown   x[2] = 3;
2509566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
251c4762a1bSJed Brown   PetscFunctionReturn(0);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown 
2549371c9d4SSatish Balay static PetscErrorCode OregoCreate(Problem p) {
255c4762a1bSJed Brown   PetscFunctionBeginUser;
256c4762a1bSJed Brown   p->destroy    = 0;
257c4762a1bSJed Brown   p->function   = &OregoFunction;
258c4762a1bSJed Brown   p->jacobian   = &OregoJacobian;
259c4762a1bSJed Brown   p->solution   = &OregoSolution;
260c4762a1bSJed Brown   p->final_time = 360;
261c4762a1bSJed Brown   p->n          = 3;
262c4762a1bSJed Brown   PetscFunctionReturn(0);
263c4762a1bSJed Brown }
264c4762a1bSJed Brown 
265c4762a1bSJed Brown /*
2660e3d61c9SBarry Smith    User-defined monitor for comparing to exact solutions when possible
267c4762a1bSJed Brown */
268c4762a1bSJed Brown typedef struct {
269c4762a1bSJed Brown   MPI_Comm comm;
270c4762a1bSJed Brown   Problem  problem;
271c4762a1bSJed Brown   Vec      x;
272c4762a1bSJed Brown } MonitorCtx;
273c4762a1bSJed Brown 
2749371c9d4SSatish Balay static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal t, Vec x, void *ctx) {
275c4762a1bSJed Brown   MonitorCtx *mon = (MonitorCtx *)ctx;
276c4762a1bSJed Brown   PetscReal   h, nrm_x, nrm_exact, nrm_diff;
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   PetscFunctionBeginUser;
279c4762a1bSJed Brown   if (!mon->problem->solution) PetscFunctionReturn(0);
2809566063dSJacob Faibussowitsch   PetscCall((*mon->problem->solution)(t, mon->x, mon->problem->data));
2819566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &nrm_x));
2829566063dSJacob Faibussowitsch   PetscCall(VecNorm(mon->x, NORM_2, &nrm_exact));
2839566063dSJacob Faibussowitsch   PetscCall(VecAYPX(mon->x, -1, x));
2849566063dSJacob Faibussowitsch   PetscCall(VecNorm(mon->x, NORM_2, &nrm_diff));
2859566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &h));
286*48a46eb9SPierre Jolivet   if (step < 0) PetscCall(PetscPrintf(mon->comm, "Interpolated final solution "));
28763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(mon->comm, "step %4" PetscInt_FMT " t=%12.8e h=% 8.2e  |x|=%9.2e  |x_e|=%9.2e  |x-x_e|=%9.2e\n", step, (double)t, (double)h, (double)nrm_x, (double)nrm_exact, (double)nrm_diff));
288c4762a1bSJed Brown   PetscFunctionReturn(0);
289c4762a1bSJed Brown }
290c4762a1bSJed Brown 
2919371c9d4SSatish Balay int main(int argc, char **argv) {
292c4762a1bSJed Brown   PetscFunctionList plist = NULL;
293c4762a1bSJed Brown   char              pname[256];
294c4762a1bSJed Brown   TS                ts;   /* nonlinear solver */
295c4762a1bSJed Brown   Vec               x, r; /* solution, residual vectors */
296c4762a1bSJed Brown   Mat               A;    /* Jacobian matrix */
297c4762a1bSJed Brown   Problem           problem;
2983d5a8a6aSBarry Smith   PetscBool         use_monitor = PETSC_FALSE;
2993d5a8a6aSBarry Smith   PetscBool         use_result  = PETSC_FALSE;
300c4762a1bSJed Brown   PetscInt          steps, nonlinits, linits, snesfails, rejects;
301c4762a1bSJed Brown   PetscReal         ftime;
302c4762a1bSJed Brown   MonitorCtx        mon;
303c4762a1bSJed Brown   PetscMPIInt       size;
304c4762a1bSJed Brown 
305c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
306c4762a1bSJed Brown      Initialize program
307c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
308327415f7SBarry Smith   PetscFunctionBeginUser;
3099566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, (char *)0, help));
3109566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3113c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /* Register the available problems */
3149566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "rober", &RoberCreate));
3159566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "ce", &CECreate));
3169566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "orego", &OregoCreate));
3179566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(pname, "ce"));
318c4762a1bSJed Brown 
319c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
320c4762a1bSJed Brown     Set runtime options
321c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
322d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Timestepping benchmark options", "");
323c4762a1bSJed Brown   {
3249566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-problem_type", "Name of problem to run", "", plist, pname, pname, sizeof(pname), NULL));
325c4762a1bSJed Brown     use_monitor = PETSC_FALSE;
3269566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-monitor_error", "Display errors relative to exact solutions", "", use_monitor, &use_monitor, NULL));
3279566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-monitor_result", "Display result", "", use_result, &use_result, NULL));
328c4762a1bSJed Brown   }
329d0609cedSBarry Smith   PetscOptionsEnd();
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   /* Create the new problem */
3329566063dSJacob Faibussowitsch   PetscCall(PetscNew(&problem));
333c4762a1bSJed Brown   problem->comm = MPI_COMM_WORLD;
334c4762a1bSJed Brown   {
335c4762a1bSJed Brown     PetscErrorCode (*pcreate)(Problem);
336c4762a1bSJed Brown 
3379566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(plist, pname, &pcreate));
3383c633725SBarry Smith     PetscCheck(pcreate, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "No problem '%s'", pname);
3399566063dSJacob Faibussowitsch     PetscCall((*pcreate)(problem));
340c4762a1bSJed Brown   }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
343c4762a1bSJed Brown     Create necessary matrix and vectors
344c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3459566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
3469566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, problem->n, problem->n, PETSC_DETERMINE, PETSC_DETERMINE));
3479566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
3489566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
349c4762a1bSJed Brown 
3509566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
3519566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
352c4762a1bSJed Brown 
353c4762a1bSJed Brown   mon.comm    = PETSC_COMM_WORLD;
354c4762a1bSJed Brown   mon.problem = problem;
3559566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &mon.x));
356c4762a1bSJed Brown 
357c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
358c4762a1bSJed Brown      Create timestepping solver context
359c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3609566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
3619566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
3629566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW)); /* Rosenbrock-W */
3639566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, problem->function, problem->data));
3649566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, problem->jacobian, problem->data));
3659566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, problem->final_time));
3669566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
3679566063dSJacob Faibussowitsch   PetscCall(TSSetMaxStepRejections(ts, 10));
3689566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSNESFailures(ts, -1)); /* unlimited */
369*48a46eb9SPierre Jolivet   if (use_monitor) PetscCall(TSMonitorSet(ts, &MonitorError, &mon, NULL));
370c4762a1bSJed Brown 
371c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
372c4762a1bSJed Brown      Set initial conditions
373c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3749566063dSJacob Faibussowitsch   PetscCall((*problem->solution)(0, x, problem->data));
3759566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .001));
3769566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
377c4762a1bSJed Brown 
378c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
379c4762a1bSJed Brown      Set runtime options
380c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3819566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
384c4762a1bSJed Brown      Solve nonlinear system
385c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3869566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
3879566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
3889566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
3899566063dSJacob Faibussowitsch   PetscCall(TSGetSNESFailures(ts, &snesfails));
3909566063dSJacob Faibussowitsch   PetscCall(TSGetStepRejections(ts, &rejects));
3919566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &nonlinits));
3929566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &linits));
3933d5a8a6aSBarry Smith   if (use_result) {
39463a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %" PetscInt_FMT " (%" PetscInt_FMT " rejected, %" PetscInt_FMT " SNES fails), ftime %g, nonlinits %" PetscInt_FMT ", linits %" PetscInt_FMT "\n", steps, rejects, snesfails, (double)ftime, nonlinits, linits));
3953d5a8a6aSBarry Smith   }
396c4762a1bSJed Brown 
397c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
398c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
399c4762a1bSJed Brown      are no longer needed.
400c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4019566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4029566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
4039566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
4049566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mon.x));
4059566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
4061baa6e33SBarry Smith   if (problem->destroy) PetscCall((*problem->destroy)(problem));
4079566063dSJacob Faibussowitsch   PetscCall(PetscFree(problem));
4089566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&plist));
409c4762a1bSJed Brown 
4109566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
411b122ec5aSJacob Faibussowitsch   return 0;
412c4762a1bSJed Brown }
413c4762a1bSJed Brown 
414c4762a1bSJed Brown /*TEST
415c4762a1bSJed Brown 
416c4762a1bSJed Brown     test:
417c4762a1bSJed Brown       requires: !complex
4183d5a8a6aSBarry Smith       args:  -monitor_result -monitor_error -ts_atol 1e-2 -ts_rtol 1e-2 -ts_exact_final_time interpolate -ts_type arkimex
419c4762a1bSJed Brown 
420c4762a1bSJed Brown     test:
421c4762a1bSJed Brown       suffix: 2
422c4762a1bSJed Brown       requires: !single !complex
4233d5a8a6aSBarry Smith       args: -monitor_result -ts_atol 1e-2 -ts_rtol 1e-2 -ts_max_time 15 -ts_type arkimex -ts_arkimex_type 2e -problem_type orego -ts_arkimex_initial_guess_extrapolate 0 -ts_adapt_time_step_increase_delay 4
424c4762a1bSJed Brown 
425c4762a1bSJed Brown     test:
426c4762a1bSJed Brown       suffix: 3
427c4762a1bSJed Brown       requires: !single !complex
4283d5a8a6aSBarry Smith       args: -monitor_result -ts_atol 1e-2 -ts_rtol 1e-2 -ts_max_time 15 -ts_type arkimex -ts_arkimex_type 2e -problem_type orego -ts_arkimex_initial_guess_extrapolate 1
4293d5a8a6aSBarry Smith 
4303d5a8a6aSBarry Smith     test:
4313d5a8a6aSBarry Smith       suffix: 4
4323d5a8a6aSBarry Smith 
4333d5a8a6aSBarry Smith     test:
4343d5a8a6aSBarry Smith       suffix: 5
4353d5a8a6aSBarry Smith       args: -snes_lag_jacobian 20 -snes_lag_jacobian_persists
436c4762a1bSJed Brown 
437c4762a1bSJed Brown TEST*/
438