xref: /petsc/src/ts/tutorials/ex8.c (revision 3d5a8a6af04c38a7eb23d2a01da917c06b399402)
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   PetscErrorCode    ierr;
34c4762a1bSJed Brown   PetscScalar       *f;
35c4762a1bSJed Brown   const PetscScalar *x,*xdot;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   PetscFunctionBeginUser;
38c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
39c4762a1bSJed Brown   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
40c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
41c4762a1bSJed Brown   f[0] = xdot[0] + 0.04*x[0] - 1e4*x[1]*x[2];
42c4762a1bSJed Brown   f[1] = xdot[1] - 0.04*x[0] + 1e4*x[1]*x[2] + 3e7*PetscSqr(x[1]);
43c4762a1bSJed Brown   f[2] = xdot[2] - 3e7*PetscSqr(x[1]);
44c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
45c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
46c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
47c4762a1bSJed Brown   PetscFunctionReturn(0);
48c4762a1bSJed Brown }
49c4762a1bSJed Brown 
50c4762a1bSJed Brown static PetscErrorCode RoberJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
51c4762a1bSJed Brown {
52c4762a1bSJed Brown   PetscErrorCode    ierr;
53c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1,2};
54c4762a1bSJed Brown   PetscScalar       J[3][3];
55c4762a1bSJed Brown   const PetscScalar *x,*xdot;
56c4762a1bSJed Brown 
57c4762a1bSJed Brown   PetscFunctionBeginUser;
58c4762a1bSJed Brown   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
59c4762a1bSJed Brown   ierr    = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
60c4762a1bSJed Brown   J[0][0] = a + 0.04;     J[0][1] = -1e4*x[2];                   J[0][2] = -1e4*x[1];
61c4762a1bSJed Brown   J[1][0] = -0.04;        J[1][1] = a + 1e4*x[2] + 3e7*2*x[1];   J[1][2] = 1e4*x[1];
62c4762a1bSJed Brown   J[2][0] = 0;            J[2][1] = -3e7*2*x[1];                 J[2][2] = a;
63c4762a1bSJed Brown   ierr    = MatSetValues(B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
64c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
65c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
66c4762a1bSJed Brown 
67c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
68c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
69c4762a1bSJed Brown   if (A != B) {
70c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
71c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
72c4762a1bSJed Brown   }
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
76c4762a1bSJed Brown static PetscErrorCode RoberSolution(PetscReal t,Vec X,void *ctx)
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   PetscErrorCode ierr;
79c4762a1bSJed Brown   PetscScalar    *x;
80c4762a1bSJed Brown 
81c4762a1bSJed Brown   PetscFunctionBeginUser;
82c4762a1bSJed Brown   if (t != 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not implemented");
83c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
84c4762a1bSJed Brown   x[0] = 1;
85c4762a1bSJed Brown   x[1] = 0;
86c4762a1bSJed Brown   x[2] = 0;
87c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
88c4762a1bSJed Brown   PetscFunctionReturn(0);
89c4762a1bSJed Brown }
90c4762a1bSJed Brown 
91c4762a1bSJed Brown static PetscErrorCode RoberCreate(Problem p)
92c4762a1bSJed Brown {
93c4762a1bSJed Brown 
94c4762a1bSJed Brown   PetscFunctionBeginUser;
95c4762a1bSJed Brown   p->destroy    = 0;
96c4762a1bSJed Brown   p->function   = &RoberFunction;
97c4762a1bSJed Brown   p->jacobian   = &RoberJacobian;
98c4762a1bSJed Brown   p->solution   = &RoberSolution;
99c4762a1bSJed Brown   p->final_time = 1e11;
100c4762a1bSJed Brown   p->n          = 3;
101c4762a1bSJed Brown   PetscFunctionReturn(0);
102c4762a1bSJed Brown }
103c4762a1bSJed Brown 
104c4762a1bSJed Brown /*
105c4762a1bSJed Brown      Stiff scalar valued problem
106c4762a1bSJed Brown */
107c4762a1bSJed Brown 
108c4762a1bSJed Brown typedef struct {
109c4762a1bSJed Brown   PetscReal lambda;
110c4762a1bSJed Brown } CECtx;
111c4762a1bSJed Brown 
112c4762a1bSJed Brown static PetscErrorCode CEDestroy(Problem p)
113c4762a1bSJed Brown {
114c4762a1bSJed Brown   PetscErrorCode ierr;
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   PetscFunctionBeginUser;
117c4762a1bSJed Brown   ierr = PetscFree(p->data);CHKERRQ(ierr);
118c4762a1bSJed Brown   PetscFunctionReturn(0);
119c4762a1bSJed Brown }
120c4762a1bSJed Brown 
121c4762a1bSJed Brown static PetscErrorCode CEFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
122c4762a1bSJed Brown {
123c4762a1bSJed Brown   PetscErrorCode    ierr;
124c4762a1bSJed Brown   PetscReal         l = ((CECtx*)ctx)->lambda;
125c4762a1bSJed Brown   PetscScalar       *f;
126c4762a1bSJed Brown   const PetscScalar *x,*xdot;
127c4762a1bSJed Brown 
128c4762a1bSJed Brown   PetscFunctionBeginUser;
129c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
130c4762a1bSJed Brown   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
131c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
132c4762a1bSJed Brown   f[0] = xdot[0] + l*(x[0] - PetscCosReal(t));
133c4762a1bSJed Brown #if 0
134c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD," f(t=%g,x=%g,xdot=%g) = %g\n",(double)t,(double)x[0],(double)xdot[0],(double)f[0]);CHKERRQ(ierr);
135c4762a1bSJed Brown #endif
136c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
137c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
138c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
139c4762a1bSJed Brown   PetscFunctionReturn(0);
140c4762a1bSJed Brown }
141c4762a1bSJed Brown 
142c4762a1bSJed Brown static PetscErrorCode CEJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
143c4762a1bSJed Brown {
144c4762a1bSJed Brown   PetscReal         l = ((CECtx*)ctx)->lambda;
145c4762a1bSJed Brown   PetscErrorCode    ierr;
146c4762a1bSJed Brown   PetscInt          rowcol[] = {0};
147c4762a1bSJed Brown   PetscScalar       J[1][1];
148c4762a1bSJed Brown   const PetscScalar *x,*xdot;
149c4762a1bSJed Brown 
150c4762a1bSJed Brown   PetscFunctionBeginUser;
151c4762a1bSJed Brown   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
152c4762a1bSJed Brown   ierr    = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
153c4762a1bSJed Brown   J[0][0] = a + l;
154c4762a1bSJed Brown   ierr    = MatSetValues(B,1,rowcol,1,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
155c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
156c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
159c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
160c4762a1bSJed Brown   if (A != B) {
161c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
162c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
163c4762a1bSJed Brown   }
164c4762a1bSJed Brown   PetscFunctionReturn(0);
165c4762a1bSJed Brown }
166c4762a1bSJed Brown 
167c4762a1bSJed Brown static PetscErrorCode CESolution(PetscReal t,Vec X,void *ctx)
168c4762a1bSJed Brown {
169c4762a1bSJed Brown   PetscReal      l = ((CECtx*)ctx)->lambda;
170c4762a1bSJed Brown   PetscErrorCode ierr;
171c4762a1bSJed Brown   PetscScalar    *x;
172c4762a1bSJed Brown 
173c4762a1bSJed Brown   PetscFunctionBeginUser;
174c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
175c4762a1bSJed Brown   x[0] = l/(l*l+1)*(l*PetscCosReal(t)+PetscSinReal(t)) - l*l/(l*l+1)*PetscExpReal(-l*t);
176c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
177c4762a1bSJed Brown   PetscFunctionReturn(0);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown static PetscErrorCode CECreate(Problem p)
181c4762a1bSJed Brown {
182c4762a1bSJed Brown   PetscErrorCode ierr;
183c4762a1bSJed Brown   CECtx          *ce;
184c4762a1bSJed Brown 
185c4762a1bSJed Brown   PetscFunctionBeginUser;
186c4762a1bSJed Brown   ierr    = PetscMalloc(sizeof(CECtx),&ce);CHKERRQ(ierr);
187c4762a1bSJed Brown   p->data = (void*)ce;
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   p->destroy    = &CEDestroy;
190c4762a1bSJed Brown   p->function   = &CEFunction;
191c4762a1bSJed Brown   p->jacobian   = &CEJacobian;
192c4762a1bSJed Brown   p->solution   = &CESolution;
193c4762a1bSJed Brown   p->final_time = 10;
194c4762a1bSJed Brown   p->n          = 1;
195c4762a1bSJed Brown   p->hasexact   = PETSC_TRUE;
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   ce->lambda = 10;
198c4762a1bSJed Brown   ierr       = PetscOptionsBegin(p->comm,NULL,"CE options","");CHKERRQ(ierr);
199c4762a1bSJed Brown   {
200c4762a1bSJed Brown     ierr = PetscOptionsReal("-problem_ce_lambda","Parameter controlling stiffness: xdot + lambda*(x - cos(t))","",ce->lambda,&ce->lambda,NULL);CHKERRQ(ierr);
201c4762a1bSJed Brown   }
202c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
203c4762a1bSJed Brown   PetscFunctionReturn(0);
204c4762a1bSJed Brown }
205c4762a1bSJed Brown 
206c4762a1bSJed Brown /*
207c4762a1bSJed Brown *  Stiff 3-variable oscillatory system from chemical reactions. problem OREGO in Hairer&Wanner
208c4762a1bSJed Brown */
209c4762a1bSJed Brown static PetscErrorCode OregoFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
210c4762a1bSJed Brown {
211c4762a1bSJed Brown   PetscErrorCode    ierr;
212c4762a1bSJed Brown   PetscScalar       *f;
213c4762a1bSJed Brown   const PetscScalar *x,*xdot;
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   PetscFunctionBeginUser;
216c4762a1bSJed Brown   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
217c4762a1bSJed Brown   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
218c4762a1bSJed Brown   ierr = VecGetArray(F,&f);CHKERRQ(ierr);
219c4762a1bSJed Brown   f[0] = xdot[0] - 77.27*(x[1] + x[0]*(1. - 8.375e-6*x[0] - x[1]));
220c4762a1bSJed Brown   f[1] = xdot[1] - 1/77.27*(x[2] - (1. + x[0])*x[1]);
221c4762a1bSJed Brown   f[2] = xdot[2] - 0.161*(x[0] - x[2]);
222c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
223c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
224c4762a1bSJed Brown   ierr = VecRestoreArray(F,&f);CHKERRQ(ierr);
225c4762a1bSJed Brown   PetscFunctionReturn(0);
226c4762a1bSJed Brown }
227c4762a1bSJed Brown 
228c4762a1bSJed Brown static PetscErrorCode OregoJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
229c4762a1bSJed Brown {
230c4762a1bSJed Brown   PetscErrorCode    ierr;
231c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1,2};
232c4762a1bSJed Brown   PetscScalar       J[3][3];
233c4762a1bSJed Brown   const PetscScalar *x,*xdot;
234c4762a1bSJed Brown 
235c4762a1bSJed Brown   PetscFunctionBeginUser;
236c4762a1bSJed Brown   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
237c4762a1bSJed Brown   ierr    = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
238c4762a1bSJed Brown   J[0][0] = a - 77.27*((1. - 8.375e-6*x[0] - x[1]) - 8.375e-6*x[0]);
239c4762a1bSJed Brown   J[0][1] = -77.27*(1. - x[0]);
240c4762a1bSJed Brown   J[0][2] = 0;
241c4762a1bSJed Brown   J[1][0] = 1./77.27*x[1];
242c4762a1bSJed Brown   J[1][1] = a + 1./77.27*(1. + x[0]);
243c4762a1bSJed Brown   J[1][2] = -1./77.27;
244c4762a1bSJed Brown   J[2][0] = -0.161;
245c4762a1bSJed Brown   J[2][1] = 0;
246c4762a1bSJed Brown   J[2][2] = a + 0.161;
247c4762a1bSJed Brown   ierr    = MatSetValues(B,3,rowcol,3,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
248c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
249c4762a1bSJed Brown   ierr    = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
250c4762a1bSJed Brown 
251c4762a1bSJed Brown   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
252c4762a1bSJed Brown   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
253c4762a1bSJed Brown   if (A != B) {
254c4762a1bSJed Brown     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
255c4762a1bSJed Brown     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
256c4762a1bSJed Brown   }
257c4762a1bSJed Brown   PetscFunctionReturn(0);
258c4762a1bSJed Brown }
259c4762a1bSJed Brown 
260c4762a1bSJed Brown static PetscErrorCode OregoSolution(PetscReal t,Vec X,void *ctx)
261c4762a1bSJed Brown {
262c4762a1bSJed Brown   PetscErrorCode ierr;
263c4762a1bSJed Brown   PetscScalar    *x;
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   PetscFunctionBeginUser;
266c4762a1bSJed Brown   if (t != 0) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"not implemented");
267c4762a1bSJed Brown   ierr = VecGetArray(X,&x);CHKERRQ(ierr);
268c4762a1bSJed Brown   x[0] = 1;
269c4762a1bSJed Brown   x[1] = 2;
270c4762a1bSJed Brown   x[2] = 3;
271c4762a1bSJed Brown   ierr = VecRestoreArray(X,&x);CHKERRQ(ierr);
272c4762a1bSJed Brown   PetscFunctionReturn(0);
273c4762a1bSJed Brown }
274c4762a1bSJed Brown 
275c4762a1bSJed Brown static PetscErrorCode OregoCreate(Problem p)
276c4762a1bSJed Brown {
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   PetscFunctionBeginUser;
279c4762a1bSJed Brown   p->destroy    = 0;
280c4762a1bSJed Brown   p->function   = &OregoFunction;
281c4762a1bSJed Brown   p->jacobian   = &OregoJacobian;
282c4762a1bSJed Brown   p->solution   = &OregoSolution;
283c4762a1bSJed Brown   p->final_time = 360;
284c4762a1bSJed Brown   p->n          = 3;
285c4762a1bSJed Brown   PetscFunctionReturn(0);
286c4762a1bSJed Brown }
287c4762a1bSJed Brown 
288c4762a1bSJed Brown 
289c4762a1bSJed Brown /*
290c4762a1bSJed Brown *  User-defined monitor for comparing to exact solutions when possible
291c4762a1bSJed Brown */
292c4762a1bSJed Brown typedef struct {
293c4762a1bSJed Brown   MPI_Comm comm;
294c4762a1bSJed Brown   Problem  problem;
295c4762a1bSJed Brown   Vec      x;
296c4762a1bSJed Brown } MonitorCtx;
297c4762a1bSJed Brown 
298c4762a1bSJed Brown static PetscErrorCode MonitorError(TS ts,PetscInt step,PetscReal t,Vec x,void *ctx)
299c4762a1bSJed Brown {
300c4762a1bSJed Brown   PetscErrorCode ierr;
301c4762a1bSJed Brown   MonitorCtx     *mon = (MonitorCtx*)ctx;
302c4762a1bSJed Brown   PetscReal      h,nrm_x,nrm_exact,nrm_diff;
303c4762a1bSJed Brown 
304c4762a1bSJed Brown   PetscFunctionBeginUser;
305c4762a1bSJed Brown   if (!mon->problem->solution) PetscFunctionReturn(0);
306c4762a1bSJed Brown   ierr = (*mon->problem->solution)(t,mon->x,mon->problem->data);CHKERRQ(ierr);
307c4762a1bSJed Brown   ierr = VecNorm(x,NORM_2,&nrm_x);CHKERRQ(ierr);
308c4762a1bSJed Brown   ierr = VecNorm(mon->x,NORM_2,&nrm_exact);CHKERRQ(ierr);
309c4762a1bSJed Brown   ierr = VecAYPX(mon->x,-1,x);CHKERRQ(ierr);
310c4762a1bSJed Brown   ierr = VecNorm(mon->x,NORM_2,&nrm_diff);CHKERRQ(ierr);
311c4762a1bSJed Brown   ierr = TSGetTimeStep(ts,&h);CHKERRQ(ierr);
312c4762a1bSJed Brown   if (step < 0) {
313c4762a1bSJed Brown     ierr = PetscPrintf(mon->comm,"Interpolated final solution ");CHKERRQ(ierr);
314c4762a1bSJed Brown   }
315c4762a1bSJed Brown   ierr = PetscPrintf(mon->comm,"step %4D 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);CHKERRQ(ierr);
316c4762a1bSJed Brown   PetscFunctionReturn(0);
317c4762a1bSJed Brown }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown 
320c4762a1bSJed Brown int main(int argc,char **argv)
321c4762a1bSJed Brown {
322c4762a1bSJed Brown   PetscFunctionList plist = NULL;
323c4762a1bSJed Brown   char              pname[256];
324c4762a1bSJed Brown   TS                ts;            /* nonlinear solver */
325c4762a1bSJed Brown   Vec               x,r;           /* solution, residual vectors */
326c4762a1bSJed Brown   Mat               A;             /* Jacobian matrix */
327c4762a1bSJed Brown   Problem           problem;
328*3d5a8a6aSBarry Smith   PetscBool         use_monitor = PETSC_FALSE;
329*3d5a8a6aSBarry Smith   PetscBool         use_result = PETSC_FALSE;
330c4762a1bSJed Brown   PetscInt          steps,nonlinits,linits,snesfails,rejects;
331c4762a1bSJed Brown   PetscReal         ftime;
332c4762a1bSJed Brown   MonitorCtx        mon;
333c4762a1bSJed Brown   PetscErrorCode    ierr;
334c4762a1bSJed Brown   PetscMPIInt       size;
335c4762a1bSJed Brown 
336c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
337c4762a1bSJed Brown      Initialize program
338c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
339c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr;
340c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
341c4762a1bSJed Brown   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs");
342c4762a1bSJed Brown 
343c4762a1bSJed Brown   /* Register the available problems */
344c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&plist,"rober",&RoberCreate);CHKERRQ(ierr);
345c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&plist,"ce",&CECreate);CHKERRQ(ierr);
346c4762a1bSJed Brown   ierr = PetscFunctionListAdd(&plist,"orego",&OregoCreate);CHKERRQ(ierr);
347c4762a1bSJed Brown   ierr = PetscStrcpy(pname,"ce");CHKERRQ(ierr);
348c4762a1bSJed Brown 
349c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
350c4762a1bSJed Brown     Set runtime options
351c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
352c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Timestepping benchmark options","");CHKERRQ(ierr);
353c4762a1bSJed Brown   {
354c4762a1bSJed Brown     ierr        = PetscOptionsFList("-problem_type","Name of problem to run","",plist,pname,pname,sizeof(pname),NULL);CHKERRQ(ierr);
355c4762a1bSJed Brown     use_monitor = PETSC_FALSE;
356c4762a1bSJed Brown     ierr        = PetscOptionsBool("-monitor_error","Display errors relative to exact solutions","",use_monitor,&use_monitor,NULL);CHKERRQ(ierr);
357*3d5a8a6aSBarry Smith     ierr        = PetscOptionsBool("-monitor_result","Display result","",use_result,&use_result,NULL);CHKERRQ(ierr);
358c4762a1bSJed Brown   }
359c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
360c4762a1bSJed Brown 
361c4762a1bSJed Brown   /* Create the new problem */
362c4762a1bSJed Brown   ierr          = PetscNew(&problem);CHKERRQ(ierr);
363c4762a1bSJed Brown   problem->comm = MPI_COMM_WORLD;
364c4762a1bSJed Brown   {
365c4762a1bSJed Brown     PetscErrorCode (*pcreate)(Problem);
366c4762a1bSJed Brown 
367c4762a1bSJed Brown     ierr = PetscFunctionListFind(plist,pname,&pcreate);CHKERRQ(ierr);
368c4762a1bSJed Brown     if (!pcreate) SETERRQ1(PETSC_COMM_SELF,1,"No problem '%s'",pname);
369c4762a1bSJed Brown     ierr = (*pcreate)(problem);CHKERRQ(ierr);
370c4762a1bSJed Brown   }
371c4762a1bSJed Brown 
372c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
373c4762a1bSJed Brown     Create necessary matrix and vectors
374c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
375c4762a1bSJed Brown   ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr);
376c4762a1bSJed Brown   ierr = MatSetSizes(A,problem->n,problem->n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
377c4762a1bSJed Brown   ierr = MatSetFromOptions(A);CHKERRQ(ierr);
378c4762a1bSJed Brown   ierr = MatSetUp(A);CHKERRQ(ierr);
379c4762a1bSJed Brown 
380c4762a1bSJed Brown   ierr = MatCreateVecs(A,&x,NULL);CHKERRQ(ierr);
381c4762a1bSJed Brown   ierr = VecDuplicate(x,&r);CHKERRQ(ierr);
382c4762a1bSJed Brown 
383c4762a1bSJed Brown   mon.comm    = PETSC_COMM_WORLD;
384c4762a1bSJed Brown   mon.problem = problem;
385c4762a1bSJed Brown   ierr        = VecDuplicate(x,&mon.x);CHKERRQ(ierr);
386c4762a1bSJed Brown 
387c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
388c4762a1bSJed Brown      Create timestepping solver context
389c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
390c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
391c4762a1bSJed Brown   ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr);
392c4762a1bSJed Brown   ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr); /* Rosenbrock-W */
393c4762a1bSJed Brown   ierr = TSSetIFunction(ts,NULL,problem->function,problem->data);CHKERRQ(ierr);
394c4762a1bSJed Brown   ierr = TSSetIJacobian(ts,A,A,problem->jacobian,problem->data);CHKERRQ(ierr);
395c4762a1bSJed Brown   ierr = TSSetMaxTime(ts,problem->final_time);CHKERRQ(ierr);
396c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
397c4762a1bSJed Brown   ierr = TSSetMaxStepRejections(ts,10);CHKERRQ(ierr);
398c4762a1bSJed Brown   ierr = TSSetMaxSNESFailures(ts,-1);CHKERRQ(ierr); /* unlimited */
399c4762a1bSJed Brown   if (use_monitor) {
400c4762a1bSJed Brown     ierr = TSMonitorSet(ts,&MonitorError,&mon,NULL);CHKERRQ(ierr);
401c4762a1bSJed Brown   }
402c4762a1bSJed Brown 
403c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
404c4762a1bSJed Brown      Set initial conditions
405c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
406c4762a1bSJed Brown   ierr = (*problem->solution)(0,x,problem->data);CHKERRQ(ierr);
407c4762a1bSJed Brown   ierr = TSSetTimeStep(ts,.001);CHKERRQ(ierr);
408c4762a1bSJed Brown   ierr = TSSetSolution(ts,x);CHKERRQ(ierr);
409c4762a1bSJed Brown 
410c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
411c4762a1bSJed Brown      Set runtime options
412c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
413c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
414c4762a1bSJed Brown 
415c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
416c4762a1bSJed Brown      Solve nonlinear system
417c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
418c4762a1bSJed Brown   ierr = TSSolve(ts,x);CHKERRQ(ierr);
419c4762a1bSJed Brown   ierr = TSGetSolveTime(ts,&ftime);CHKERRQ(ierr);
420c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&steps);CHKERRQ(ierr);
421c4762a1bSJed Brown   ierr = TSGetSNESFailures(ts,&snesfails);CHKERRQ(ierr);
422c4762a1bSJed Brown   ierr = TSGetStepRejections(ts,&rejects);CHKERRQ(ierr);
423c4762a1bSJed Brown   ierr = TSGetSNESIterations(ts,&nonlinits);CHKERRQ(ierr);
424c4762a1bSJed Brown   ierr = TSGetKSPIterations(ts,&linits);CHKERRQ(ierr);
425*3d5a8a6aSBarry Smith   if (use_result) {
426c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD,"steps %D (%D rejected, %D SNES fails), ftime %g, nonlinits %D, linits %D\n",steps,rejects,snesfails,(double)ftime,nonlinits,linits);CHKERRQ(ierr);
427*3d5a8a6aSBarry Smith   }
428c4762a1bSJed Brown 
429c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
430c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
431c4762a1bSJed Brown      are no longer needed.
432c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
433c4762a1bSJed Brown   ierr = MatDestroy(&A);CHKERRQ(ierr);
434c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
435c4762a1bSJed Brown   ierr = VecDestroy(&r);CHKERRQ(ierr);
436c4762a1bSJed Brown   ierr = VecDestroy(&mon.x);CHKERRQ(ierr);
437c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
438c4762a1bSJed Brown   if (problem->destroy) {
439c4762a1bSJed Brown     ierr = (*problem->destroy)(problem);CHKERRQ(ierr);
440c4762a1bSJed Brown   }
441c4762a1bSJed Brown   ierr = PetscFree(problem);CHKERRQ(ierr);
442c4762a1bSJed Brown   ierr = PetscFunctionListDestroy(&plist);CHKERRQ(ierr);
443c4762a1bSJed Brown 
444c4762a1bSJed Brown   ierr = PetscFinalize();
445c4762a1bSJed Brown   return ierr;
446c4762a1bSJed Brown }
447c4762a1bSJed Brown 
448c4762a1bSJed Brown /*TEST
449c4762a1bSJed Brown 
450c4762a1bSJed Brown     test:
451c4762a1bSJed Brown       requires: !complex
452*3d5a8a6aSBarry Smith       args:  -monitor_result -monitor_error -ts_atol 1e-2 -ts_rtol 1e-2 -ts_exact_final_time interpolate -ts_type arkimex
453c4762a1bSJed Brown 
454c4762a1bSJed Brown     test:
455c4762a1bSJed Brown       suffix: 2
456c4762a1bSJed Brown       requires: !single !complex
457*3d5a8a6aSBarry 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
458c4762a1bSJed Brown 
459c4762a1bSJed Brown     test:
460c4762a1bSJed Brown       suffix: 3
461c4762a1bSJed Brown       requires: !single !complex
462*3d5a8a6aSBarry 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
463*3d5a8a6aSBarry Smith 
464*3d5a8a6aSBarry Smith     test:
465*3d5a8a6aSBarry Smith       suffix: 4
466*3d5a8a6aSBarry Smith 
467*3d5a8a6aSBarry Smith     test:
468*3d5a8a6aSBarry Smith       suffix: 5
469*3d5a8a6aSBarry Smith       args: -snes_lag_jacobian 20 -snes_lag_jacobian_persists
470c4762a1bSJed Brown 
471c4762a1bSJed Brown TEST*/
472