xref: /petsc/src/ts/tutorials/ex8.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 */
31c4762a1bSJed Brown static PetscErrorCode RoberFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
32c4762a1bSJed Brown {
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));
46c4762a1bSJed Brown   PetscFunctionReturn(0);
47c4762a1bSJed Brown }
48c4762a1bSJed Brown 
49c4762a1bSJed Brown static PetscErrorCode RoberJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
50c4762a1bSJed Brown {
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));
58c4762a1bSJed Brown   J[0][0] = a + 0.04;     J[0][1] = -1e4*x[2];                   J[0][2] = -1e4*x[1];
59c4762a1bSJed Brown   J[1][0] = -0.04;        J[1][1] = a + 1e4*x[2] + 3e7*2*x[1];   J[1][2] = 1e4*x[1];
60c4762a1bSJed Brown   J[2][0] = 0;            J[2][1] = -3e7*2*x[1];                 J[2][2] = a;
619566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES));
629566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
639566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot,&xdot));
64c4762a1bSJed Brown 
659566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
669566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
67c4762a1bSJed Brown   if (A != B) {
689566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
699566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
70c4762a1bSJed Brown   }
71c4762a1bSJed Brown   PetscFunctionReturn(0);
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
74c4762a1bSJed Brown static PetscErrorCode RoberSolution(PetscReal t,Vec X,void *ctx)
75c4762a1bSJed Brown {
76c4762a1bSJed Brown   PetscScalar    *x;
77c4762a1bSJed Brown 
78c4762a1bSJed Brown   PetscFunctionBeginUser;
793c633725SBarry Smith   PetscCheck(t == 0,PETSC_COMM_WORLD,PETSC_ERR_SUP,"not implemented");
809566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
81c4762a1bSJed Brown   x[0] = 1;
82c4762a1bSJed Brown   x[1] = 0;
83c4762a1bSJed Brown   x[2] = 0;
849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
85c4762a1bSJed Brown   PetscFunctionReturn(0);
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown static PetscErrorCode RoberCreate(Problem p)
89c4762a1bSJed Brown {
90c4762a1bSJed Brown   PetscFunctionBeginUser;
91c4762a1bSJed Brown   p->destroy    = 0;
92c4762a1bSJed Brown   p->function   = &RoberFunction;
93c4762a1bSJed Brown   p->jacobian   = &RoberJacobian;
94c4762a1bSJed Brown   p->solution   = &RoberSolution;
95c4762a1bSJed Brown   p->final_time = 1e11;
96c4762a1bSJed Brown   p->n          = 3;
97c4762a1bSJed Brown   PetscFunctionReturn(0);
98c4762a1bSJed Brown }
99c4762a1bSJed Brown 
100c4762a1bSJed Brown /*
101c4762a1bSJed Brown      Stiff scalar valued problem
102c4762a1bSJed Brown */
103c4762a1bSJed Brown 
104c4762a1bSJed Brown typedef struct {
105c4762a1bSJed Brown   PetscReal lambda;
106c4762a1bSJed Brown } CECtx;
107c4762a1bSJed Brown 
108c4762a1bSJed Brown static PetscErrorCode CEDestroy(Problem p)
109c4762a1bSJed Brown {
110c4762a1bSJed Brown   PetscFunctionBeginUser;
1119566063dSJacob Faibussowitsch   PetscCall(PetscFree(p->data));
112c4762a1bSJed Brown   PetscFunctionReturn(0);
113c4762a1bSJed Brown }
114c4762a1bSJed Brown 
115c4762a1bSJed Brown static PetscErrorCode CEFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
116c4762a1bSJed Brown {
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 
135c4762a1bSJed Brown static PetscErrorCode CEJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
136c4762a1bSJed Brown {
137c4762a1bSJed Brown   PetscReal         l = ((CECtx*)ctx)->lambda;
138c4762a1bSJed Brown   PetscInt          rowcol[] = {0};
139c4762a1bSJed Brown   PetscScalar       J[1][1];
140c4762a1bSJed Brown   const PetscScalar *x,*xdot;
141c4762a1bSJed Brown 
142c4762a1bSJed Brown   PetscFunctionBeginUser;
1439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
1449566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot,&xdot));
145c4762a1bSJed Brown   J[0][0] = a + l;
1469566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,1,rowcol,1,rowcol,&J[0][0],INSERT_VALUES));
1479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
1489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot,&xdot));
149c4762a1bSJed Brown 
1509566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1519566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
152c4762a1bSJed Brown   if (A != B) {
1539566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1549566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
155c4762a1bSJed Brown   }
156c4762a1bSJed Brown   PetscFunctionReturn(0);
157c4762a1bSJed Brown }
158c4762a1bSJed Brown 
159c4762a1bSJed Brown static PetscErrorCode CESolution(PetscReal t,Vec X,void *ctx)
160c4762a1bSJed Brown {
161c4762a1bSJed Brown   PetscReal      l = ((CECtx*)ctx)->lambda;
162c4762a1bSJed Brown   PetscScalar    *x;
163c4762a1bSJed Brown 
164c4762a1bSJed Brown   PetscFunctionBeginUser;
1659566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
166c4762a1bSJed Brown   x[0] = l/(l*l+1)*(l*PetscCosReal(t)+PetscSinReal(t)) - l*l/(l*l+1)*PetscExpReal(-l*t);
1679566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
168c4762a1bSJed Brown   PetscFunctionReturn(0);
169c4762a1bSJed Brown }
170c4762a1bSJed Brown 
171c4762a1bSJed Brown static PetscErrorCode CECreate(Problem p)
172c4762a1bSJed Brown {
173c4762a1bSJed Brown   CECtx          *ce;
174c4762a1bSJed Brown 
175c4762a1bSJed Brown   PetscFunctionBeginUser;
1769566063dSJacob Faibussowitsch   PetscCall(PetscMalloc(sizeof(CECtx),&ce));
177c4762a1bSJed Brown   p->data = (void*)ce;
178c4762a1bSJed Brown 
179c4762a1bSJed Brown   p->destroy    = &CEDestroy;
180c4762a1bSJed Brown   p->function   = &CEFunction;
181c4762a1bSJed Brown   p->jacobian   = &CEJacobian;
182c4762a1bSJed Brown   p->solution   = &CESolution;
183c4762a1bSJed Brown   p->final_time = 10;
184c4762a1bSJed Brown   p->n          = 1;
185c4762a1bSJed Brown   p->hasexact   = PETSC_TRUE;
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   ce->lambda = 10;
188d0609cedSBarry Smith   PetscOptionsBegin(p->comm,NULL,"CE options","");
189c4762a1bSJed Brown   {
1909566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-problem_ce_lambda","Parameter controlling stiffness: xdot + lambda*(x - cos(t))","",ce->lambda,&ce->lambda,NULL));
191c4762a1bSJed Brown   }
192d0609cedSBarry Smith   PetscOptionsEnd();
193c4762a1bSJed Brown   PetscFunctionReturn(0);
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown /*
1970e3d61c9SBarry Smith    Stiff 3-variable oscillatory system from chemical reactions. problem OREGO in Hairer&Wanner
198c4762a1bSJed Brown */
199c4762a1bSJed Brown static PetscErrorCode OregoFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
200c4762a1bSJed Brown {
201c4762a1bSJed Brown   PetscScalar       *f;
202c4762a1bSJed Brown   const PetscScalar *x,*xdot;
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   PetscFunctionBeginUser;
2059566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
2069566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot,&xdot));
2079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
208c4762a1bSJed Brown   f[0] = xdot[0] - 77.27*(x[1] + x[0]*(1. - 8.375e-6*x[0] - x[1]));
209c4762a1bSJed Brown   f[1] = xdot[1] - 1/77.27*(x[2] - (1. + x[0])*x[1]);
210c4762a1bSJed Brown   f[2] = xdot[2] - 0.161*(x[0] - x[2]);
2119566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
2129566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot,&xdot));
2139566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
214c4762a1bSJed Brown   PetscFunctionReturn(0);
215c4762a1bSJed Brown }
216c4762a1bSJed Brown 
217c4762a1bSJed Brown static PetscErrorCode OregoJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
218c4762a1bSJed Brown {
219c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1,2};
220c4762a1bSJed Brown   PetscScalar       J[3][3];
221c4762a1bSJed Brown   const PetscScalar *x,*xdot;
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   PetscFunctionBeginUser;
2249566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
2259566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot,&xdot));
226c4762a1bSJed Brown   J[0][0] = a - 77.27*((1. - 8.375e-6*x[0] - x[1]) - 8.375e-6*x[0]);
227c4762a1bSJed Brown   J[0][1] = -77.27*(1. - x[0]);
228c4762a1bSJed Brown   J[0][2] = 0;
229c4762a1bSJed Brown   J[1][0] = 1./77.27*x[1];
230c4762a1bSJed Brown   J[1][1] = a + 1./77.27*(1. + x[0]);
231c4762a1bSJed Brown   J[1][2] = -1./77.27;
232c4762a1bSJed Brown   J[2][0] = -0.161;
233c4762a1bSJed Brown   J[2][1] = 0;
234c4762a1bSJed Brown   J[2][2] = a + 0.161;
2359566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES));
2369566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
2379566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot,&xdot));
238c4762a1bSJed Brown 
2399566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
2409566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
241c4762a1bSJed Brown   if (A != B) {
2429566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
2439566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
244c4762a1bSJed Brown   }
245c4762a1bSJed Brown   PetscFunctionReturn(0);
246c4762a1bSJed Brown }
247c4762a1bSJed Brown 
248c4762a1bSJed Brown static PetscErrorCode OregoSolution(PetscReal t,Vec X,void *ctx)
249c4762a1bSJed Brown {
250c4762a1bSJed Brown   PetscScalar    *x;
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   PetscFunctionBeginUser;
2533c633725SBarry Smith   PetscCheck(t == 0,PETSC_COMM_WORLD,PETSC_ERR_SUP,"not implemented");
2549566063dSJacob Faibussowitsch   PetscCall(VecGetArray(X,&x));
255c4762a1bSJed Brown   x[0] = 1;
256c4762a1bSJed Brown   x[1] = 2;
257c4762a1bSJed Brown   x[2] = 3;
2589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(X,&x));
259c4762a1bSJed Brown   PetscFunctionReturn(0);
260c4762a1bSJed Brown }
261c4762a1bSJed Brown 
262c4762a1bSJed Brown static PetscErrorCode OregoCreate(Problem p)
263c4762a1bSJed Brown {
264c4762a1bSJed Brown   PetscFunctionBeginUser;
265c4762a1bSJed Brown   p->destroy    = 0;
266c4762a1bSJed Brown   p->function   = &OregoFunction;
267c4762a1bSJed Brown   p->jacobian   = &OregoJacobian;
268c4762a1bSJed Brown   p->solution   = &OregoSolution;
269c4762a1bSJed Brown   p->final_time = 360;
270c4762a1bSJed Brown   p->n          = 3;
271c4762a1bSJed Brown   PetscFunctionReturn(0);
272c4762a1bSJed Brown }
273c4762a1bSJed Brown 
274c4762a1bSJed Brown /*
2750e3d61c9SBarry Smith    User-defined monitor for comparing to exact solutions when possible
276c4762a1bSJed Brown */
277c4762a1bSJed Brown typedef struct {
278c4762a1bSJed Brown   MPI_Comm comm;
279c4762a1bSJed Brown   Problem  problem;
280c4762a1bSJed Brown   Vec      x;
281c4762a1bSJed Brown } MonitorCtx;
282c4762a1bSJed Brown 
283c4762a1bSJed Brown static PetscErrorCode MonitorError(TS ts,PetscInt step,PetscReal t,Vec x,void *ctx)
284c4762a1bSJed Brown {
285c4762a1bSJed Brown   MonitorCtx     *mon = (MonitorCtx*)ctx;
286c4762a1bSJed Brown   PetscReal      h,nrm_x,nrm_exact,nrm_diff;
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   PetscFunctionBeginUser;
289c4762a1bSJed Brown   if (!mon->problem->solution) PetscFunctionReturn(0);
2909566063dSJacob Faibussowitsch   PetscCall((*mon->problem->solution)(t,mon->x,mon->problem->data));
2919566063dSJacob Faibussowitsch   PetscCall(VecNorm(x,NORM_2,&nrm_x));
2929566063dSJacob Faibussowitsch   PetscCall(VecNorm(mon->x,NORM_2,&nrm_exact));
2939566063dSJacob Faibussowitsch   PetscCall(VecAYPX(mon->x,-1,x));
2949566063dSJacob Faibussowitsch   PetscCall(VecNorm(mon->x,NORM_2,&nrm_diff));
2959566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts,&h));
296c4762a1bSJed Brown   if (step < 0) {
2979566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(mon->comm,"Interpolated final solution "));
298c4762a1bSJed Brown   }
29963a3b9bcSJacob 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));
300c4762a1bSJed Brown   PetscFunctionReturn(0);
301c4762a1bSJed Brown }
302c4762a1bSJed Brown 
303c4762a1bSJed Brown int main(int argc,char **argv)
304c4762a1bSJed Brown {
305c4762a1bSJed Brown   PetscFunctionList plist = NULL;
306c4762a1bSJed Brown   char              pname[256];
307c4762a1bSJed Brown   TS                ts;            /* nonlinear solver */
308c4762a1bSJed Brown   Vec               x,r;           /* solution, residual vectors */
309c4762a1bSJed Brown   Mat               A;             /* Jacobian matrix */
310c4762a1bSJed Brown   Problem           problem;
3113d5a8a6aSBarry Smith   PetscBool         use_monitor = PETSC_FALSE;
3123d5a8a6aSBarry Smith   PetscBool         use_result = PETSC_FALSE;
313c4762a1bSJed Brown   PetscInt          steps,nonlinits,linits,snesfails,rejects;
314c4762a1bSJed Brown   PetscReal         ftime;
315c4762a1bSJed Brown   MonitorCtx        mon;
316c4762a1bSJed Brown   PetscMPIInt       size;
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
319c4762a1bSJed Brown      Initialize program
320c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
321*327415f7SBarry Smith   PetscFunctionBeginUser;
3229566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,(char*)0,help));
3239566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
3243c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   /* Register the available problems */
3279566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist,"rober",&RoberCreate));
3289566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist,"ce",&CECreate));
3299566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListAdd(&plist,"orego",&OregoCreate));
3309566063dSJacob Faibussowitsch   PetscCall(PetscStrcpy(pname,"ce"));
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
333c4762a1bSJed Brown     Set runtime options
334c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
335d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Timestepping benchmark options","");
336c4762a1bSJed Brown   {
3379566063dSJacob Faibussowitsch     PetscCall(PetscOptionsFList("-problem_type","Name of problem to run","",plist,pname,pname,sizeof(pname),NULL));
338c4762a1bSJed Brown     use_monitor = PETSC_FALSE;
3399566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-monitor_error","Display errors relative to exact solutions","",use_monitor,&use_monitor,NULL));
3409566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-monitor_result","Display result","",use_result,&use_result,NULL));
341c4762a1bSJed Brown   }
342d0609cedSBarry Smith   PetscOptionsEnd();
343c4762a1bSJed Brown 
344c4762a1bSJed Brown   /* Create the new problem */
3459566063dSJacob Faibussowitsch   PetscCall(PetscNew(&problem));
346c4762a1bSJed Brown   problem->comm = MPI_COMM_WORLD;
347c4762a1bSJed Brown   {
348c4762a1bSJed Brown     PetscErrorCode (*pcreate)(Problem);
349c4762a1bSJed Brown 
3509566063dSJacob Faibussowitsch     PetscCall(PetscFunctionListFind(plist,pname,&pcreate));
3513c633725SBarry Smith     PetscCheck(pcreate,PETSC_COMM_SELF,PETSC_ERR_ARG_UNKNOWN_TYPE,"No problem '%s'",pname);
3529566063dSJacob Faibussowitsch     PetscCall((*pcreate)(problem));
353c4762a1bSJed Brown   }
354c4762a1bSJed Brown 
355c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
356c4762a1bSJed Brown     Create necessary matrix and vectors
357c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3589566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&A));
3599566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A,problem->n,problem->n,PETSC_DETERMINE,PETSC_DETERMINE));
3609566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
3619566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
362c4762a1bSJed Brown 
3639566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A,&x,NULL));
3649566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&r));
365c4762a1bSJed Brown 
366c4762a1bSJed Brown   mon.comm    = PETSC_COMM_WORLD;
367c4762a1bSJed Brown   mon.problem = problem;
3689566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(x,&mon.x));
369c4762a1bSJed Brown 
370c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
371c4762a1bSJed Brown      Create timestepping solver context
372c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3739566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
3749566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts,TS_NONLINEAR));
3759566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSROSW)); /* Rosenbrock-W */
3769566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,problem->function,problem->data));
3779566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,A,A,problem->jacobian,problem->data));
3789566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,problem->final_time));
3799566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
3809566063dSJacob Faibussowitsch   PetscCall(TSSetMaxStepRejections(ts,10));
3819566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSNESFailures(ts,-1)); /* unlimited */
382c4762a1bSJed Brown   if (use_monitor) {
3839566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,&MonitorError,&mon,NULL));
384c4762a1bSJed Brown   }
385c4762a1bSJed Brown 
386c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
387c4762a1bSJed Brown      Set initial conditions
388c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3899566063dSJacob Faibussowitsch   PetscCall((*problem->solution)(0,x,problem->data));
3909566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,.001));
3919566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts,x));
392c4762a1bSJed Brown 
393c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
394c4762a1bSJed Brown      Set runtime options
395c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
3969566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
397c4762a1bSJed Brown 
398c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
399c4762a1bSJed Brown      Solve nonlinear system
400c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4019566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,x));
4029566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&ftime));
4039566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&steps));
4049566063dSJacob Faibussowitsch   PetscCall(TSGetSNESFailures(ts,&snesfails));
4059566063dSJacob Faibussowitsch   PetscCall(TSGetStepRejections(ts,&rejects));
4069566063dSJacob Faibussowitsch   PetscCall(TSGetSNESIterations(ts,&nonlinits));
4079566063dSJacob Faibussowitsch   PetscCall(TSGetKSPIterations(ts,&linits));
4083d5a8a6aSBarry Smith   if (use_result) {
40963a3b9bcSJacob 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));
4103d5a8a6aSBarry Smith   }
411c4762a1bSJed Brown 
412c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
413c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
414c4762a1bSJed Brown      are no longer needed.
415c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
4169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
4179566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
4189566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&r));
4199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&mon.x));
4209566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
4211baa6e33SBarry Smith   if (problem->destroy) PetscCall((*problem->destroy)(problem));
4229566063dSJacob Faibussowitsch   PetscCall(PetscFree(problem));
4239566063dSJacob Faibussowitsch   PetscCall(PetscFunctionListDestroy(&plist));
424c4762a1bSJed Brown 
4259566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
426b122ec5aSJacob Faibussowitsch   return 0;
427c4762a1bSJed Brown }
428c4762a1bSJed Brown 
429c4762a1bSJed Brown /*TEST
430c4762a1bSJed Brown 
431c4762a1bSJed Brown     test:
432c4762a1bSJed Brown       requires: !complex
4333d5a8a6aSBarry Smith       args:  -monitor_result -monitor_error -ts_atol 1e-2 -ts_rtol 1e-2 -ts_exact_final_time interpolate -ts_type arkimex
434c4762a1bSJed Brown 
435c4762a1bSJed Brown     test:
436c4762a1bSJed Brown       suffix: 2
437c4762a1bSJed Brown       requires: !single !complex
4383d5a8a6aSBarry 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
439c4762a1bSJed Brown 
440c4762a1bSJed Brown     test:
441c4762a1bSJed Brown       suffix: 3
442c4762a1bSJed Brown       requires: !single !complex
4433d5a8a6aSBarry 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
4443d5a8a6aSBarry Smith 
4453d5a8a6aSBarry Smith     test:
4463d5a8a6aSBarry Smith       suffix: 4
4473d5a8a6aSBarry Smith 
4483d5a8a6aSBarry Smith     test:
4493d5a8a6aSBarry Smith       suffix: 5
4503d5a8a6aSBarry Smith       args: -snes_lag_jacobian 20 -snes_lag_jacobian_persists
451c4762a1bSJed Brown 
452c4762a1bSJed Brown TEST*/
453