xref: /petsc/src/ts/tutorials/ex8.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "Nonlinear DAE benchmark problems.\n";
2c4762a1bSJed Brown 
3c4762a1bSJed Brown /*
4c4762a1bSJed Brown    Include "petscts.h" so that we can use TS solvers.  Note that this
5c4762a1bSJed Brown    file automatically includes:
6c4762a1bSJed Brown      petscsys.h       - base PETSc routines   petscvec.h - vectors
7c4762a1bSJed Brown      petscmat.h - matrices
8c4762a1bSJed Brown      petscis.h     - index sets            petscksp.h - Krylov subspace methods
9c4762a1bSJed Brown      petscviewer.h - viewers               petscpc.h  - preconditioners
10c4762a1bSJed Brown      petscksp.h   - linear solvers
11c4762a1bSJed Brown */
12c4762a1bSJed Brown #include <petscts.h>
13c4762a1bSJed Brown 
14c4762a1bSJed Brown typedef struct _Problem *Problem;
15c4762a1bSJed Brown struct _Problem {
16c4762a1bSJed Brown   PetscErrorCode (*destroy)(Problem);
178434afd1SBarry Smith   TSIFunctionFn *function;
188434afd1SBarry Smith   TSIJacobianFn *jacobian;
19c4762a1bSJed Brown   PetscErrorCode (*solution)(PetscReal, Vec, void *);
20c4762a1bSJed Brown   MPI_Comm  comm;
21c4762a1bSJed Brown   PetscReal final_time;
22c4762a1bSJed Brown   PetscInt  n;
23c4762a1bSJed Brown   PetscBool hasexact;
24c4762a1bSJed Brown   void     *data;
25c4762a1bSJed Brown };
26c4762a1bSJed Brown 
27c4762a1bSJed Brown /*
28c4762a1bSJed Brown       Stiff 3-variable system from chemical reactions, due to Robertson (1966), problem ROBER in Hairer&Wanner, ODE 2, 1996
29efa39862SBarry Smith       https://archimede.uniba.it/~testset/report/rober.pdf
30c4762a1bSJed Brown */
RoberFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)31*2a8381b2SBarry Smith static PetscErrorCode RoberFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx)
32d71ae5a4SJacob Faibussowitsch {
33c4762a1bSJed Brown   PetscScalar       *f;
34c4762a1bSJed Brown   const PetscScalar *x, *xdot;
35c4762a1bSJed Brown 
36c4762a1bSJed Brown   PetscFunctionBeginUser;
379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
40c4762a1bSJed Brown   f[0] = xdot[0] + 0.04 * x[0] - 1e4 * x[1] * x[2];
41c4762a1bSJed Brown   f[1] = xdot[1] - 0.04 * x[0] + 1e4 * x[1] * x[2] + 3e7 * PetscSqr(x[1]);
42c4762a1bSJed Brown   f[2] = xdot[2] - 3e7 * PetscSqr(x[1]);
439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
449566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
459566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
RoberJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,PetscCtx ctx)49*2a8381b2SBarry Smith static PetscErrorCode RoberJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
50d71ae5a4SJacob Faibussowitsch {
51c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1, 2};
52c4762a1bSJed Brown   PetscScalar        J[3][3];
53c4762a1bSJed Brown   const PetscScalar *x, *xdot;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown   PetscFunctionBeginUser;
569566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
579566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
589371c9d4SSatish Balay   J[0][0] = a + 0.04;
599371c9d4SSatish Balay   J[0][1] = -1e4 * x[2];
609371c9d4SSatish Balay   J[0][2] = -1e4 * x[1];
619371c9d4SSatish Balay   J[1][0] = -0.04;
629371c9d4SSatish Balay   J[1][1] = a + 1e4 * x[2] + 3e7 * 2 * x[1];
639371c9d4SSatish Balay   J[1][2] = 1e4 * x[1];
649371c9d4SSatish Balay   J[2][0] = 0;
659371c9d4SSatish Balay   J[2][1] = -3e7 * 2 * x[1];
669371c9d4SSatish Balay   J[2][2] = a;
679566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 3, rowcol, 3, rowcol, &J[0][0], INSERT_VALUES));
689566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
70c4762a1bSJed Brown 
719566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
729566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
73c4762a1bSJed Brown   if (A != B) {
749566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
759566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
76c4762a1bSJed Brown   }
773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
78c4762a1bSJed Brown }
79c4762a1bSJed Brown 
RoberSolution(PetscReal t,Vec X,PetscCtx ctx)80*2a8381b2SBarry Smith static PetscErrorCode RoberSolution(PetscReal t, Vec X, PetscCtx ctx)
81d71ae5a4SJacob Faibussowitsch {
82c4762a1bSJed Brown   PetscScalar *x;
83c4762a1bSJed Brown 
84c4762a1bSJed Brown   PetscFunctionBeginUser;
853c633725SBarry Smith   PetscCheck(t == 0, PETSC_COMM_WORLD, PETSC_ERR_SUP, "not implemented");
869566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
87c4762a1bSJed Brown   x[0] = 1;
88c4762a1bSJed Brown   x[1] = 0;
89c4762a1bSJed Brown   x[2] = 0;
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
913ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
92c4762a1bSJed Brown }
93c4762a1bSJed Brown 
RoberCreate(Problem p)94d71ae5a4SJacob Faibussowitsch static PetscErrorCode RoberCreate(Problem p)
95d71ae5a4SJacob Faibussowitsch {
96c4762a1bSJed Brown   PetscFunctionBeginUser;
97c4762a1bSJed Brown   p->destroy    = 0;
98c4762a1bSJed Brown   p->function   = &RoberFunction;
99c4762a1bSJed Brown   p->jacobian   = &RoberJacobian;
100c4762a1bSJed Brown   p->solution   = &RoberSolution;
101c4762a1bSJed Brown   p->final_time = 1e11;
102c4762a1bSJed Brown   p->n          = 3;
1033ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
104c4762a1bSJed Brown }
105c4762a1bSJed Brown 
106c4762a1bSJed Brown /*
107c4762a1bSJed Brown      Stiff scalar valued problem
108c4762a1bSJed Brown */
109c4762a1bSJed Brown 
110c4762a1bSJed Brown typedef struct {
111c4762a1bSJed Brown   PetscReal lambda;
112c4762a1bSJed Brown } CECtx;
113c4762a1bSJed Brown 
CEDestroy(Problem p)114d71ae5a4SJacob Faibussowitsch static PetscErrorCode CEDestroy(Problem p)
115d71ae5a4SJacob Faibussowitsch {
116c4762a1bSJed Brown   PetscFunctionBeginUser;
1179566063dSJacob Faibussowitsch   PetscCall(PetscFree(p->data));
1183ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
119c4762a1bSJed Brown }
120c4762a1bSJed Brown 
CEFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)121*2a8381b2SBarry Smith static PetscErrorCode CEFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx)
122d71ae5a4SJacob Faibussowitsch {
123c4762a1bSJed Brown   PetscReal          l = ((CECtx *)ctx)->lambda;
124c4762a1bSJed Brown   PetscScalar       *f;
125c4762a1bSJed Brown   const PetscScalar *x, *xdot;
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   PetscFunctionBeginUser;
1289566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1299566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
1309566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
131c4762a1bSJed Brown   f[0] = xdot[0] + l * (x[0] - PetscCosReal(t));
132c4762a1bSJed Brown #if 0
1339566063dSJacob 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]));
134c4762a1bSJed Brown #endif
1359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
1379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
1383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
CEJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,PetscCtx ctx)141*2a8381b2SBarry Smith static PetscErrorCode CEJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
142d71ae5a4SJacob Faibussowitsch {
143c4762a1bSJed Brown   PetscReal          l        = ((CECtx *)ctx)->lambda;
144c4762a1bSJed Brown   PetscInt           rowcol[] = {0};
145c4762a1bSJed Brown   PetscScalar        J[1][1];
146c4762a1bSJed Brown   const PetscScalar *x, *xdot;
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   PetscFunctionBeginUser;
1499566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
151c4762a1bSJed Brown   J[0][0] = a + l;
1529566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 1, rowcol, 1, rowcol, &J[0][0], INSERT_VALUES));
1539566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1549566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
155c4762a1bSJed Brown 
1569566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1579566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
158c4762a1bSJed Brown   if (A != B) {
1599566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1609566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
161c4762a1bSJed Brown   }
1623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
163c4762a1bSJed Brown }
164c4762a1bSJed Brown 
CESolution(PetscReal t,Vec X,PetscCtx ctx)165*2a8381b2SBarry Smith static PetscErrorCode CESolution(PetscReal t, Vec X, PetscCtx ctx)
166d71ae5a4SJacob Faibussowitsch {
167c4762a1bSJed Brown   PetscReal    l = ((CECtx *)ctx)->lambda;
168c4762a1bSJed Brown   PetscScalar *x;
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   PetscFunctionBeginUser;
1719566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
172c4762a1bSJed Brown   x[0] = l / (l * l + 1) * (l * PetscCosReal(t) + PetscSinReal(t)) - l * l / (l * l + 1) * PetscExpReal(-l * t);
1739566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
1743ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
175c4762a1bSJed Brown }
176c4762a1bSJed Brown 
CECreate(Problem p)177d71ae5a4SJacob Faibussowitsch static PetscErrorCode CECreate(Problem p)
178d71ae5a4SJacob Faibussowitsch {
179c4762a1bSJed Brown   CECtx *ce;
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   PetscFunctionBeginUser;
1829566063dSJacob Faibussowitsch   PetscCall(PetscMalloc(sizeof(CECtx), &ce));
183c4762a1bSJed Brown   p->data = (void *)ce;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   p->destroy    = &CEDestroy;
186c4762a1bSJed Brown   p->function   = &CEFunction;
187c4762a1bSJed Brown   p->jacobian   = &CEJacobian;
188c4762a1bSJed Brown   p->solution   = &CESolution;
189c4762a1bSJed Brown   p->final_time = 10;
190c4762a1bSJed Brown   p->n          = 1;
191c4762a1bSJed Brown   p->hasexact   = PETSC_TRUE;
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   ce->lambda = 10;
194d0609cedSBarry Smith   PetscOptionsBegin(p->comm, NULL, "CE options", "");
195d71ae5a4SJacob Faibussowitsch   {
196d71ae5a4SJacob Faibussowitsch     PetscCall(PetscOptionsReal("-problem_ce_lambda", "Parameter controlling stiffness: xdot + lambda*(x - cos(t))", "", ce->lambda, &ce->lambda, NULL));
197d71ae5a4SJacob Faibussowitsch   }
198d0609cedSBarry Smith   PetscOptionsEnd();
1993ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
200c4762a1bSJed Brown }
201c4762a1bSJed Brown 
202c4762a1bSJed Brown /*
2030e3d61c9SBarry Smith    Stiff 3-variable oscillatory system from chemical reactions. problem OREGO in Hairer&Wanner
204c4762a1bSJed Brown */
OregoFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)205*2a8381b2SBarry Smith static PetscErrorCode OregoFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx)
206d71ae5a4SJacob Faibussowitsch {
207c4762a1bSJed Brown   PetscScalar       *f;
208c4762a1bSJed Brown   const PetscScalar *x, *xdot;
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   PetscFunctionBeginUser;
2119566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
2139566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
214c4762a1bSJed Brown   f[0] = xdot[0] - 77.27 * (x[1] + x[0] * (1. - 8.375e-6 * x[0] - x[1]));
215c4762a1bSJed Brown   f[1] = xdot[1] - 1 / 77.27 * (x[2] - (1. + x[0]) * x[1]);
216c4762a1bSJed Brown   f[2] = xdot[2] - 0.161 * (x[0] - x[2]);
2179566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2189566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
2199566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
2203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
221c4762a1bSJed Brown }
222c4762a1bSJed Brown 
OregoJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,PetscCtx ctx)223*2a8381b2SBarry Smith static PetscErrorCode OregoJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
224d71ae5a4SJacob Faibussowitsch {
225c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1, 2};
226c4762a1bSJed Brown   PetscScalar        J[3][3];
227c4762a1bSJed Brown   const PetscScalar *x, *xdot;
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   PetscFunctionBeginUser;
2309566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
2319566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
232c4762a1bSJed Brown   J[0][0] = a - 77.27 * ((1. - 8.375e-6 * x[0] - x[1]) - 8.375e-6 * x[0]);
233c4762a1bSJed Brown   J[0][1] = -77.27 * (1. - x[0]);
234c4762a1bSJed Brown   J[0][2] = 0;
235c4762a1bSJed Brown   J[1][0] = 1. / 77.27 * x[1];
236c4762a1bSJed Brown   J[1][1] = a + 1. / 77.27 * (1. + x[0]);
237c4762a1bSJed Brown   J[1][2] = -1. / 77.27;
238c4762a1bSJed Brown   J[2][0] = -0.161;
239c4762a1bSJed Brown   J[2][1] = 0;
240c4762a1bSJed Brown   J[2][2] = a + 0.161;
2419566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 3, rowcol, 3, rowcol, &J[0][0], INSERT_VALUES));
2429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
2439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
244c4762a1bSJed Brown 
2459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
2469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
247c4762a1bSJed Brown   if (A != B) {
2489566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
2499566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
250c4762a1bSJed Brown   }
2513ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
252c4762a1bSJed Brown }
253c4762a1bSJed Brown 
OregoSolution(PetscReal t,Vec X,PetscCtx ctx)254*2a8381b2SBarry Smith static PetscErrorCode OregoSolution(PetscReal t, Vec X, PetscCtx ctx)
255d71ae5a4SJacob Faibussowitsch {
256c4762a1bSJed Brown   PetscScalar *x;
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   PetscFunctionBeginUser;
2593c633725SBarry Smith   PetscCheck(t == 0, PETSC_COMM_WORLD, PETSC_ERR_SUP, "not implemented");
2609566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X, &x));
261c4762a1bSJed Brown   x[0] = 1;
262c4762a1bSJed Brown   x[1] = 2;
263c4762a1bSJed Brown   x[2] = 3;
2649566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X, &x));
2653ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
266c4762a1bSJed Brown }
267c4762a1bSJed Brown 
OregoCreate(Problem p)268d71ae5a4SJacob Faibussowitsch static PetscErrorCode OregoCreate(Problem p)
269d71ae5a4SJacob Faibussowitsch {
270c4762a1bSJed Brown   PetscFunctionBeginUser;
271c4762a1bSJed Brown   p->destroy    = 0;
272c4762a1bSJed Brown   p->function   = &OregoFunction;
273c4762a1bSJed Brown   p->jacobian   = &OregoJacobian;
274c4762a1bSJed Brown   p->solution   = &OregoSolution;
275c4762a1bSJed Brown   p->final_time = 360;
276c4762a1bSJed Brown   p->n          = 3;
2773ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
278c4762a1bSJed Brown }
279c4762a1bSJed Brown 
280c4762a1bSJed Brown /*
2810e3d61c9SBarry Smith    User-defined monitor for comparing to exact solutions when possible
282c4762a1bSJed Brown */
283c4762a1bSJed Brown typedef struct {
284c4762a1bSJed Brown   MPI_Comm comm;
285c4762a1bSJed Brown   Problem  problem;
286c4762a1bSJed Brown   Vec      x;
287c4762a1bSJed Brown } MonitorCtx;
288c4762a1bSJed Brown 
MonitorError(TS ts,PetscInt step,PetscReal t,Vec x,PetscCtx ctx)289*2a8381b2SBarry Smith static PetscErrorCode MonitorError(TS ts, PetscInt step, PetscReal t, Vec x, PetscCtx ctx)
290d71ae5a4SJacob Faibussowitsch {
291c4762a1bSJed Brown   MonitorCtx *mon = (MonitorCtx *)ctx;
292c4762a1bSJed Brown   PetscReal   h, nrm_x, nrm_exact, nrm_diff;
293c4762a1bSJed Brown 
294c4762a1bSJed Brown   PetscFunctionBeginUser;
2953ba16761SJacob Faibussowitsch   if (!mon->problem->solution) PetscFunctionReturn(PETSC_SUCCESS);
2969566063dSJacob Faibussowitsch   PetscCall((*mon->problem->solution)(t, mon->x, mon->problem->data));
2979566063dSJacob Faibussowitsch   PetscCall(VecNorm(x, NORM_2, &nrm_x));
2989566063dSJacob Faibussowitsch   PetscCall(VecNorm(mon->x, NORM_2, &nrm_exact));
2999566063dSJacob Faibussowitsch   PetscCall(VecAYPX(mon->x, -1, x));
3009566063dSJacob Faibussowitsch   PetscCall(VecNorm(mon->x, NORM_2, &nrm_diff));
3019566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &h));
30248a46eb9SPierre Jolivet   if (step < 0) PetscCall(PetscPrintf(mon->comm, "Interpolated final solution "));
30363a3b9bcSJacob 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));
3043ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
305c4762a1bSJed Brown }
306c4762a1bSJed Brown 
main(int argc,char ** argv)307d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
308d71ae5a4SJacob Faibussowitsch {
309c4762a1bSJed Brown   PetscFunctionList plist = NULL;
310c4762a1bSJed Brown   char              pname[256];
311c4762a1bSJed Brown   TS                ts;   /* nonlinear solver */
312c4762a1bSJed Brown   Vec               x, r; /* solution, residual vectors */
313c4762a1bSJed Brown   Mat               A;    /* Jacobian matrix */
314c4762a1bSJed Brown   Problem           problem;
3153d5a8a6aSBarry Smith   PetscBool         use_monitor = PETSC_FALSE;
3163d5a8a6aSBarry Smith   PetscBool         use_result  = PETSC_FALSE;
317c4762a1bSJed Brown   PetscInt          steps, nonlinits, linits, snesfails, rejects;
318c4762a1bSJed Brown   PetscReal         ftime;
319c4762a1bSJed Brown   MonitorCtx        mon;
320c4762a1bSJed Brown   PetscMPIInt       size;
321c4762a1bSJed Brown 
322c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
323c4762a1bSJed Brown      Initialize program
324c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
325327415f7SBarry Smith   PetscFunctionBeginUser;
326c8025a54SPierre Jolivet   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
3279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
3283c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs");
329c4762a1bSJed Brown 
330c4762a1bSJed Brown   /* Register the available problems */
3319566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "rober", &RoberCreate));
3329566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "ce", &CECreate));
3339566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist, "orego", &OregoCreate));
334c6a7a370SJeremy L Thompson   PetscCall(PetscStrncpy(pname, "ce", sizeof(pname)));
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
337c4762a1bSJed Brown     Set runtime options
338c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
339d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Timestepping benchmark options", "");
340c4762a1bSJed Brown   {
3419566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-problem_type", "Name of problem to run", "", plist, pname, pname, sizeof(pname), NULL));
342c4762a1bSJed Brown     use_monitor = PETSC_FALSE;
3439566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-monitor_error", "Display errors relative to exact solutions", "", use_monitor, &use_monitor, NULL));
3449566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-monitor_result", "Display result", "", use_result, &use_result, NULL));
345c4762a1bSJed Brown   }
346d0609cedSBarry Smith   PetscOptionsEnd();
347c4762a1bSJed Brown 
348c4762a1bSJed Brown   /* Create the new problem */
3499566063dSJacob Faibussowitsch   PetscCall(PetscNew(&problem));
350c4762a1bSJed Brown   problem->comm = MPI_COMM_WORLD;
351c4762a1bSJed Brown   {
352c4762a1bSJed Brown     PetscErrorCode (*pcreate)(Problem);
353c4762a1bSJed Brown 
3549566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(plist, pname, &pcreate));
3553c633725SBarry Smith     PetscCheck(pcreate, PETSC_COMM_SELF, PETSC_ERR_ARG_UNKNOWN_TYPE, "No problem '%s'", pname);
3569566063dSJacob Faibussowitsch     PetscCall((*pcreate)(problem));
357c4762a1bSJed Brown   }
358c4762a1bSJed Brown 
359c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
360c4762a1bSJed Brown     Create necessary matrix and vectors
361c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3629566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
3639566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, problem->n, problem->n, PETSC_DETERMINE, PETSC_DETERMINE));
3649566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
3659566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
366c4762a1bSJed Brown 
3679566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
3689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &r));
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   mon.comm    = PETSC_COMM_WORLD;
371c4762a1bSJed Brown   mon.problem = problem;
3729566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x, &mon.x));
373c4762a1bSJed Brown 
374c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
375c4762a1bSJed Brown      Create timestepping solver context
376c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3779566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
3789566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
3799566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSROSW)); /* Rosenbrock-W */
3809566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, problem->function, problem->data));
3819566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, problem->jacobian, problem->data));
3829566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, problem->final_time));
3839566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
3849566063dSJacob Faibussowitsch   PetscCall(TSSetMaxStepRejections(ts, 10));
3859566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSNESFailures(ts, -1)); /* unlimited */
38648a46eb9SPierre Jolivet   if (use_monitor) PetscCall(TSMonitorSet(ts, &MonitorError, &mon, NULL));
387c4762a1bSJed Brown 
388c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
389c4762a1bSJed Brown      Set initial conditions
390c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3919566063dSJacob Faibussowitsch   PetscCall((*problem->solution)(0, x, problem->data));
3929566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .001));
3939566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
394c4762a1bSJed Brown 
395c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
396c4762a1bSJed Brown      Set runtime options
397c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3989566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
399c4762a1bSJed Brown 
400c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
401c4762a1bSJed Brown      Solve nonlinear system
402c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4039566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
4049566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
4059566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
4069566063dSJacob Faibussowitsch   PetscCall(TSGetSNESFailures(ts, &snesfails));
4079566063dSJacob Faibussowitsch   PetscCall(TSGetStepRejections(ts, &rejects));
4089566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts, &nonlinits));
4099566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts, &linits));
4103d5a8a6aSBarry Smith   if (use_result) {
41163a3b9bcSJacob 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));
4123d5a8a6aSBarry Smith   }
413c4762a1bSJed Brown 
414c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
415c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
416c4762a1bSJed Brown      are no longer needed.
417c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
4209566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
4219566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mon.x));
4229566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
4231baa6e33SBarry Smith   if (problem->destroy) PetscCall((*problem->destroy)(problem));
4249566063dSJacob Faibussowitsch   PetscCall(PetscFree(problem));
4259566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&plist));
426c4762a1bSJed Brown 
4279566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
428b122ec5aSJacob Faibussowitsch   return 0;
429c4762a1bSJed Brown }
430c4762a1bSJed Brown 
431c4762a1bSJed Brown /*TEST
432c4762a1bSJed Brown 
433c4762a1bSJed Brown     test:
434c4762a1bSJed Brown       requires: !complex
4353d5a8a6aSBarry Smith       args: -monitor_result -monitor_error -ts_atol 1e-2 -ts_rtol 1e-2 -ts_exact_final_time interpolate -ts_type arkimex
436c4762a1bSJed Brown 
437c4762a1bSJed Brown     test:
438c4762a1bSJed Brown       suffix: 2
439c4762a1bSJed Brown       requires: !single !complex
4403d5a8a6aSBarry 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
441c4762a1bSJed Brown 
442c4762a1bSJed Brown     test:
443c4762a1bSJed Brown       suffix: 3
444c4762a1bSJed Brown       requires: !single !complex
4453d5a8a6aSBarry 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
4463d5a8a6aSBarry Smith 
4473d5a8a6aSBarry Smith     test:
4483d5a8a6aSBarry Smith       suffix: 4
4493886731fSPierre Jolivet       output_file: output/empty.out
4503d5a8a6aSBarry Smith 
4513d5a8a6aSBarry Smith     test:
4523d5a8a6aSBarry Smith       suffix: 5
4533d5a8a6aSBarry Smith       args: -snes_lag_jacobian 20 -snes_lag_jacobian_persists
4543886731fSPierre Jolivet       output_file: output/empty.out
455c4762a1bSJed Brown 
456c4762a1bSJed Brown TEST*/
457