xref: /petsc/src/ts/tutorials/ex42.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "Meinhard't activator-inhibitor model to test TS domain error feature.\n";
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown /*
4*c4762a1bSJed Brown    The activator-inhibitor on a line is described by the PDE:
5*c4762a1bSJed Brown 
6*c4762a1bSJed Brown    da/dt = \alpha a^2 / (1 + \beta h) + \rho_a - \mu_a a + D_a d^2 a/ dx^2
7*c4762a1bSJed Brown    dh/dt = \alpha a^2 + \rho_h - \mu_h h + D_h d^2 h/ dx^2
8*c4762a1bSJed Brown 
9*c4762a1bSJed Brown    The PDE part will be solve by finite-difference on the line of cells.
10*c4762a1bSJed Brown  */
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown #include <petscts.h>
13*c4762a1bSJed Brown 
14*c4762a1bSJed Brown typedef struct {
15*c4762a1bSJed Brown   PetscInt  nb_cells;
16*c4762a1bSJed Brown   PetscReal alpha;
17*c4762a1bSJed Brown   PetscReal beta;
18*c4762a1bSJed Brown   PetscReal rho_a;
19*c4762a1bSJed Brown   PetscReal rho_h;
20*c4762a1bSJed Brown   PetscReal mu_a;
21*c4762a1bSJed Brown   PetscReal mu_h;
22*c4762a1bSJed Brown   PetscReal D_a;
23*c4762a1bSJed Brown   PetscReal D_h;
24*c4762a1bSJed Brown } AppCtx;
25*c4762a1bSJed Brown 
26*c4762a1bSJed Brown PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec DXDT, void* ptr)
27*c4762a1bSJed Brown {
28*c4762a1bSJed Brown   AppCtx*           user = (AppCtx*)ptr;
29*c4762a1bSJed Brown   PetscInt          nb_cells, i;
30*c4762a1bSJed Brown   PetscReal         alpha, beta;
31*c4762a1bSJed Brown   PetscReal         rho_a, mu_a, D_a;
32*c4762a1bSJed Brown   PetscReal         rho_h, mu_h, D_h;
33*c4762a1bSJed Brown   PetscReal         a, h, da, dh, d2a, d2h;
34*c4762a1bSJed Brown   PetscErrorCode    ierr;
35*c4762a1bSJed Brown   PetscScalar       *dxdt;
36*c4762a1bSJed Brown   const PetscScalar *x;
37*c4762a1bSJed Brown 
38*c4762a1bSJed Brown   PetscFunctionBegin;
39*c4762a1bSJed Brown   nb_cells = user->nb_cells;
40*c4762a1bSJed Brown   alpha    = user->alpha;
41*c4762a1bSJed Brown   beta     = user->beta;
42*c4762a1bSJed Brown   rho_a    = user->rho_a;
43*c4762a1bSJed Brown   mu_a     = user->mu_a;
44*c4762a1bSJed Brown   D_a      = user->D_a;
45*c4762a1bSJed Brown   rho_h    = user->rho_h;
46*c4762a1bSJed Brown   mu_h     = user->mu_h;
47*c4762a1bSJed Brown   D_h      = user->D_h;
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
50*c4762a1bSJed Brown   ierr = VecGetArray(DXDT, &dxdt);CHKERRQ(ierr);
51*c4762a1bSJed Brown 
52*c4762a1bSJed Brown   for(i = 0 ; i < nb_cells ; i++) {
53*c4762a1bSJed Brown     a = x[2*i];
54*c4762a1bSJed Brown     h = x[2*i+1];
55*c4762a1bSJed Brown     // Reaction:
56*c4762a1bSJed Brown     da = alpha * a*a / (1. + beta * h) + rho_a - mu_a * a;
57*c4762a1bSJed Brown     dh = alpha * a*a + rho_h - mu_h*h;
58*c4762a1bSJed Brown     // Diffusion:
59*c4762a1bSJed Brown     d2a = d2h = 0.;
60*c4762a1bSJed Brown     if(i > 0) {
61*c4762a1bSJed Brown       d2a += (x[2*(i-1)] - a);
62*c4762a1bSJed Brown       d2h += (x[2*(i-1)+1] - h);
63*c4762a1bSJed Brown     }
64*c4762a1bSJed Brown     if(i < nb_cells-1) {
65*c4762a1bSJed Brown       d2a += (x[2*(i+1)] - a);
66*c4762a1bSJed Brown       d2h += (x[2*(i+1)+1] - h);
67*c4762a1bSJed Brown     }
68*c4762a1bSJed Brown     dxdt[2*i] = da + D_a*d2a;
69*c4762a1bSJed Brown     dxdt[2*i+1] = dh + D_h*d2h;
70*c4762a1bSJed Brown   }
71*c4762a1bSJed Brown   ierr = VecRestoreArray(DXDT, &dxdt);CHKERRQ(ierr);
72*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
73*c4762a1bSJed Brown   PetscFunctionReturn(0);
74*c4762a1bSJed Brown }
75*c4762a1bSJed Brown 
76*c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat J, Mat B, void *ptr)
77*c4762a1bSJed Brown {
78*c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
79*c4762a1bSJed Brown   PetscInt          nb_cells, i, idx;
80*c4762a1bSJed Brown   PetscReal         alpha, beta;
81*c4762a1bSJed Brown   PetscReal         mu_a, D_a;
82*c4762a1bSJed Brown   PetscReal         mu_h, D_h;
83*c4762a1bSJed Brown   PetscReal         a, h;
84*c4762a1bSJed Brown   const PetscScalar *x;
85*c4762a1bSJed Brown   PetscScalar       va[4], vh[4];
86*c4762a1bSJed Brown   PetscInt          ca[4], ch[4], rowa, rowh;
87*c4762a1bSJed Brown   PetscErrorCode    ierr;
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown   PetscFunctionBegin;
90*c4762a1bSJed Brown   nb_cells = user->nb_cells;
91*c4762a1bSJed Brown   alpha    = user->alpha;
92*c4762a1bSJed Brown   beta     = user->beta;
93*c4762a1bSJed Brown   mu_a     = user->mu_a;
94*c4762a1bSJed Brown   D_a      = user->D_a;
95*c4762a1bSJed Brown   mu_h     = user->mu_h;
96*c4762a1bSJed Brown   D_h      = user->D_h;
97*c4762a1bSJed Brown 
98*c4762a1bSJed Brown   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
99*c4762a1bSJed Brown   for(i = 0; i < nb_cells ; ++i) {
100*c4762a1bSJed Brown     rowa = 2*i;
101*c4762a1bSJed Brown     rowh = 2*i+1;
102*c4762a1bSJed Brown     a = x[2*i];
103*c4762a1bSJed Brown     h = x[2*i+1];
104*c4762a1bSJed Brown     ca[0] = ch[1] = 2*i;
105*c4762a1bSJed Brown     va[0] = 2*alpha*a / (1.+beta*h) - mu_a;
106*c4762a1bSJed Brown     vh[1] = 2*alpha*a;
107*c4762a1bSJed Brown     ca[1] = ch[0] = 2*i+1;
108*c4762a1bSJed Brown     va[1] = -beta*alpha*a*a / ((1.+beta*h)*(1.+beta*h));
109*c4762a1bSJed Brown     vh[0] = -mu_h;
110*c4762a1bSJed Brown     idx = 2;
111*c4762a1bSJed Brown     if(i > 0) {
112*c4762a1bSJed Brown       ca[idx] = 2*(i-1);
113*c4762a1bSJed Brown       ch[idx] = 2*(i-1)+1;
114*c4762a1bSJed Brown       va[idx] = D_a;
115*c4762a1bSJed Brown       vh[idx] = D_h;
116*c4762a1bSJed Brown       va[0] -= D_a;
117*c4762a1bSJed Brown       vh[0] -= D_h;
118*c4762a1bSJed Brown       idx++;
119*c4762a1bSJed Brown     }
120*c4762a1bSJed Brown     if(i < nb_cells-1) {
121*c4762a1bSJed Brown       ca[idx] = 2*(i+1);
122*c4762a1bSJed Brown       ch[idx] = 2*(i+1)+1;
123*c4762a1bSJed Brown       va[idx] = D_a;
124*c4762a1bSJed Brown       vh[idx] = D_h;
125*c4762a1bSJed Brown       va[0] -= D_a;
126*c4762a1bSJed Brown       vh[0] -= D_h;
127*c4762a1bSJed Brown       idx++;
128*c4762a1bSJed Brown     }
129*c4762a1bSJed Brown     ierr = MatSetValues(B, 1, &rowa, idx, ca, va, INSERT_VALUES);CHKERRQ(ierr);
130*c4762a1bSJed Brown     ierr = MatSetValues(B, 1, &rowh, idx, ch, vh, INSERT_VALUES);CHKERRQ(ierr);
131*c4762a1bSJed Brown   }
132*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
133*c4762a1bSJed Brown   ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
134*c4762a1bSJed Brown   ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
135*c4762a1bSJed Brown   if (J != B) {
136*c4762a1bSJed Brown     ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
137*c4762a1bSJed Brown     ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
138*c4762a1bSJed Brown   }
139*c4762a1bSJed Brown   PetscFunctionReturn(0);
140*c4762a1bSJed Brown }
141*c4762a1bSJed Brown 
142*c4762a1bSJed Brown PetscErrorCode DomainErrorFunction(TS ts, PetscReal t, Vec Y, PetscBool *accept)
143*c4762a1bSJed Brown {
144*c4762a1bSJed Brown   AppCtx            *user;
145*c4762a1bSJed Brown   PetscReal         dt;
146*c4762a1bSJed Brown   PetscErrorCode    ierr;
147*c4762a1bSJed Brown   const PetscScalar *x;
148*c4762a1bSJed Brown   PetscInt          nb_cells, i;
149*c4762a1bSJed Brown 
150*c4762a1bSJed Brown   ierr = TSGetApplicationContext(ts, &user);CHKERRQ(ierr);
151*c4762a1bSJed Brown   nb_cells = user->nb_cells;
152*c4762a1bSJed Brown   ierr = VecGetArrayRead(Y, &x);CHKERRQ(ierr);
153*c4762a1bSJed Brown   for(i = 0 ; i < 2*nb_cells ; ++i) {
154*c4762a1bSJed Brown     if(PetscRealPart(x[i]) < 0) {
155*c4762a1bSJed Brown       ierr = TSGetTimeStep(ts, &dt);CHKERRQ(ierr);
156*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_WORLD, " ** Domain Error at time %g\n", (double)t);CHKERRQ(ierr);
157*c4762a1bSJed Brown       *accept = PETSC_FALSE;
158*c4762a1bSJed Brown       break;
159*c4762a1bSJed Brown     }
160*c4762a1bSJed Brown   }
161*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(Y, &x);CHKERRQ(ierr);
162*c4762a1bSJed Brown   PetscFunctionReturn(0);
163*c4762a1bSJed Brown }
164*c4762a1bSJed Brown 
165*c4762a1bSJed Brown PetscErrorCode FormInitialState(Vec X, AppCtx* user)
166*c4762a1bSJed Brown {
167*c4762a1bSJed Brown   PetscErrorCode ierr;
168*c4762a1bSJed Brown   PetscRandom    R;
169*c4762a1bSJed Brown 
170*c4762a1bSJed Brown   PetscFunctionBegin;
171*c4762a1bSJed Brown   ierr = PetscRandomCreate(PETSC_COMM_WORLD, &R);CHKERRQ(ierr);
172*c4762a1bSJed Brown   ierr = PetscRandomSetFromOptions(R);CHKERRQ(ierr);
173*c4762a1bSJed Brown   ierr = PetscRandomSetInterval(R, 0., 10.);CHKERRQ(ierr);
174*c4762a1bSJed Brown 
175*c4762a1bSJed Brown   /*
176*c4762a1bSJed Brown    * Initialize the state vector
177*c4762a1bSJed Brown    */
178*c4762a1bSJed Brown   ierr = VecSetRandom(X, R);CHKERRQ(ierr);
179*c4762a1bSJed Brown   ierr = PetscRandomDestroy(&R);CHKERRQ(ierr);
180*c4762a1bSJed Brown   PetscFunctionReturn(0);
181*c4762a1bSJed Brown }
182*c4762a1bSJed Brown 
183*c4762a1bSJed Brown PetscErrorCode PrintSolution(Vec X, AppCtx *user)
184*c4762a1bSJed Brown {
185*c4762a1bSJed Brown   PetscErrorCode    ierr;
186*c4762a1bSJed Brown   const PetscScalar *x;
187*c4762a1bSJed Brown   PetscInt          i;
188*c4762a1bSJed Brown   PetscInt          nb_cells = user->nb_cells;
189*c4762a1bSJed Brown 
190*c4762a1bSJed Brown   PetscFunctionBegin;
191*c4762a1bSJed Brown   ierr = VecGetArrayRead(X, &x);CHKERRQ(ierr);
192*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "Activator,Inhibitor\n");CHKERRQ(ierr);
193*c4762a1bSJed Brown   for(i = 0 ; i < nb_cells ; i++) {
194*c4762a1bSJed Brown     ierr = PetscPrintf(PETSC_COMM_WORLD, "%5.6e,%5.6e\n", (double)x[2*i], (double)x[2*i+1]);CHKERRQ(ierr);
195*c4762a1bSJed Brown   }
196*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(X, &x);CHKERRQ(ierr);
197*c4762a1bSJed Brown   PetscFunctionReturn(0);
198*c4762a1bSJed Brown }
199*c4762a1bSJed Brown 
200*c4762a1bSJed Brown int main(int argc, char **argv)
201*c4762a1bSJed Brown {
202*c4762a1bSJed Brown   TS             ts;       /* time-stepping context */
203*c4762a1bSJed Brown   Vec            x;       /* State vector */
204*c4762a1bSJed Brown   Mat            J; /* Jacobian matrix */
205*c4762a1bSJed Brown   AppCtx         user; /* user-defined context */
206*c4762a1bSJed Brown   PetscErrorCode ierr;
207*c4762a1bSJed Brown   PetscReal      ftime;
208*c4762a1bSJed Brown   PetscInt       its;
209*c4762a1bSJed Brown   PetscMPIInt    size;
210*c4762a1bSJed Brown 
211*c4762a1bSJed Brown   ierr = PetscInitialize(&argc, &argv, NULL, help);if (ierr) return ierr;
212*c4762a1bSJed Brown   ierr = MPI_Comm_size(PETSC_COMM_WORLD, &size);CHKERRQ(ierr);
213*c4762a1bSJed Brown   if(size != 1) SETERRQ(PETSC_COMM_WORLD, PETSC_ERR_SUP, "This is a uniprocessor example only");
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown   /*
216*c4762a1bSJed Brown    * Allow user to set the grid dimensions and the equations parameters
217*c4762a1bSJed Brown    */
218*c4762a1bSJed Brown 
219*c4762a1bSJed Brown   user.nb_cells = 50;
220*c4762a1bSJed Brown   user.alpha = 10.;
221*c4762a1bSJed Brown   user.beta = 1.;
222*c4762a1bSJed Brown   user.rho_a = 1.;
223*c4762a1bSJed Brown   user.rho_h = 2.;
224*c4762a1bSJed Brown   user.mu_a = 2.;
225*c4762a1bSJed Brown   user.mu_h = 3.;
226*c4762a1bSJed Brown   user.D_a = 0.;
227*c4762a1bSJed Brown   user.D_h = 30.;
228*c4762a1bSJed Brown 
229*c4762a1bSJed Brown   ierr = PetscOptionsBegin(PETSC_COMM_WORLD, "", "Problem settings", "PROBLEM");CHKERRQ(ierr);
230*c4762a1bSJed Brown   ierr = PetscOptionsInt("-nb_cells", "Number of cells", "ex42.c",user.nb_cells, &user.nb_cells,NULL);CHKERRQ(ierr);
231*c4762a1bSJed Brown   ierr = PetscOptionsReal("-alpha", "Autocatalysis factor", "ex42.c",user.alpha, &user.alpha,NULL);CHKERRQ(ierr);
232*c4762a1bSJed Brown   ierr = PetscOptionsReal("-beta", "Inhibition factor", "ex42.c",user.beta, &user.beta,NULL);CHKERRQ(ierr);
233*c4762a1bSJed Brown   ierr = PetscOptionsReal("-rho_a", "Default production of the activator", "ex42.c",user.rho_a, &user.rho_a,NULL);CHKERRQ(ierr);
234*c4762a1bSJed Brown   ierr = PetscOptionsReal("-mu_a", "Degradation rate of the activator", "ex42.c",user.mu_a, &user.mu_a,NULL);CHKERRQ(ierr);
235*c4762a1bSJed Brown   ierr = PetscOptionsReal("-D_a", "Diffusion rate of the activator", "ex42.c",user.D_a, &user.D_a,NULL);CHKERRQ(ierr);
236*c4762a1bSJed Brown   ierr = PetscOptionsReal("-rho_h", "Default production of the inhibitor", "ex42.c",user.rho_h, &user.rho_h,NULL);CHKERRQ(ierr);
237*c4762a1bSJed Brown   ierr = PetscOptionsReal("-mu_h", "Degradation rate of the inhibitor", "ex42.c",user.mu_h, &user.mu_h,NULL);CHKERRQ(ierr);
238*c4762a1bSJed Brown   ierr = PetscOptionsReal("-D_h", "Diffusion rate of the inhibitor", "ex42.c",user.D_h, &user.D_h,NULL);CHKERRQ(ierr);
239*c4762a1bSJed Brown   ierr = PetscOptionsEnd();CHKERRQ(ierr);
240*c4762a1bSJed Brown 
241*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "nb_cells: %D\n", user.nb_cells);CHKERRQ(ierr);
242*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "alpha: %5.5g\n", (double)user.alpha);CHKERRQ(ierr);
243*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "beta:  %5.5g\n", (double)user.beta);CHKERRQ(ierr);
244*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "rho_a: %5.5g\n", (double)user.rho_a);CHKERRQ(ierr);
245*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "mu_a:  %5.5g\n", (double)user.mu_a);CHKERRQ(ierr);
246*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "D_a:   %5.5g\n", (double)user.D_a);CHKERRQ(ierr);
247*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "rho_h: %5.5g\n", (double)user.rho_h);CHKERRQ(ierr);
248*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "mu_h:  %5.5g\n", (double)user.mu_h);CHKERRQ(ierr);
249*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "D_h:   %5.5g\n", (double)user.D_h);CHKERRQ(ierr);
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown   /*
252*c4762a1bSJed Brown    * Create vector to hold the solution
253*c4762a1bSJed Brown    */
254*c4762a1bSJed Brown   ierr = VecCreateSeq(PETSC_COMM_WORLD, 2*user.nb_cells, &x);CHKERRQ(ierr);
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown   /*
257*c4762a1bSJed Brown    * Create time-stepper context
258*c4762a1bSJed Brown    */
259*c4762a1bSJed Brown   ierr = TSCreate(PETSC_COMM_WORLD, &ts);CHKERRQ(ierr);
260*c4762a1bSJed Brown   ierr = TSSetProblemType(ts, TS_NONLINEAR);CHKERRQ(ierr);
261*c4762a1bSJed Brown 
262*c4762a1bSJed Brown   /*
263*c4762a1bSJed Brown    * Tell the time-stepper context where to compute the solution
264*c4762a1bSJed Brown    */
265*c4762a1bSJed Brown   ierr = TSSetSolution(ts, x);CHKERRQ(ierr);
266*c4762a1bSJed Brown 
267*c4762a1bSJed Brown   /*
268*c4762a1bSJed Brown    * Allocate the jacobian matrix
269*c4762a1bSJed Brown    */
270*c4762a1bSJed Brown   ierr = MatCreateSeqAIJ(PETSC_COMM_WORLD, 2*user.nb_cells, 2*user.nb_cells, 4, 0, &J);CHKERRQ(ierr);
271*c4762a1bSJed Brown 
272*c4762a1bSJed Brown   /*
273*c4762a1bSJed Brown    * Provide the call-back for the non-linear function we are evaluating.
274*c4762a1bSJed Brown    */
275*c4762a1bSJed Brown   ierr = TSSetRHSFunction(ts, NULL, RHSFunction, &user);CHKERRQ(ierr);
276*c4762a1bSJed Brown 
277*c4762a1bSJed Brown   /*
278*c4762a1bSJed Brown    * Set the Jacobian matrix and the function user to compute Jacobians
279*c4762a1bSJed Brown    */
280*c4762a1bSJed Brown   ierr = TSSetRHSJacobian(ts, J, J, RHSJacobian, &user);CHKERRQ(ierr);
281*c4762a1bSJed Brown 
282*c4762a1bSJed Brown   /*
283*c4762a1bSJed Brown    * Set the function checking the domain
284*c4762a1bSJed Brown    */
285*c4762a1bSJed Brown   ierr = TSSetFunctionDomainError(ts, &DomainErrorFunction);CHKERRQ(ierr);
286*c4762a1bSJed Brown 
287*c4762a1bSJed Brown   /*
288*c4762a1bSJed Brown    * Initialize the problem with random values
289*c4762a1bSJed Brown    */
290*c4762a1bSJed Brown   ierr = FormInitialState(x, &user);CHKERRQ(ierr);
291*c4762a1bSJed Brown 
292*c4762a1bSJed Brown   /*
293*c4762a1bSJed Brown    * Read the solver type from options
294*c4762a1bSJed Brown    */
295*c4762a1bSJed Brown   ierr = TSSetType(ts, TSPSEUDO);CHKERRQ(ierr);
296*c4762a1bSJed Brown 
297*c4762a1bSJed Brown   /*
298*c4762a1bSJed Brown    * Set a large number of timesteps and final duration time to insure
299*c4762a1bSJed Brown    * convergenge to steady state
300*c4762a1bSJed Brown    */
301*c4762a1bSJed Brown   ierr = TSSetMaxSteps(ts, 2147483647);CHKERRQ(ierr);
302*c4762a1bSJed Brown   ierr = TSSetMaxTime(ts, 1.e12);CHKERRQ(ierr);
303*c4762a1bSJed Brown   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);CHKERRQ(ierr);
304*c4762a1bSJed Brown 
305*c4762a1bSJed Brown   /*
306*c4762a1bSJed Brown    * Set a larger number of potential errors
307*c4762a1bSJed Brown    */
308*c4762a1bSJed Brown   ierr = TSSetMaxStepRejections(ts, 50);CHKERRQ(ierr);
309*c4762a1bSJed Brown 
310*c4762a1bSJed Brown   /*
311*c4762a1bSJed Brown    * Also start with a very small dt
312*c4762a1bSJed Brown    */
313*c4762a1bSJed Brown   ierr = TSSetTimeStep(ts, 0.05);CHKERRQ(ierr);
314*c4762a1bSJed Brown 
315*c4762a1bSJed Brown   /*
316*c4762a1bSJed Brown    * Set a larger time step increment
317*c4762a1bSJed Brown    */
318*c4762a1bSJed Brown   ierr = TSPseudoSetTimeStepIncrement(ts, 1.5);CHKERRQ(ierr);
319*c4762a1bSJed Brown 
320*c4762a1bSJed Brown   /*
321*c4762a1bSJed Brown    * Let the user personalise TS
322*c4762a1bSJed Brown    */
323*c4762a1bSJed Brown   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
324*c4762a1bSJed Brown 
325*c4762a1bSJed Brown   /*
326*c4762a1bSJed Brown    * Set the context for the time stepper
327*c4762a1bSJed Brown    */
328*c4762a1bSJed Brown   ierr = TSSetApplicationContext(ts, &user);CHKERRQ(ierr);
329*c4762a1bSJed Brown 
330*c4762a1bSJed Brown   /*
331*c4762a1bSJed Brown    * Setup the time stepper, ready for evaluation
332*c4762a1bSJed Brown    */
333*c4762a1bSJed Brown   ierr = TSSetUp(ts);CHKERRQ(ierr);
334*c4762a1bSJed Brown 
335*c4762a1bSJed Brown   /*
336*c4762a1bSJed Brown    * Perform the solve.
337*c4762a1bSJed Brown    */
338*c4762a1bSJed Brown   ierr = TSSolve(ts, x);CHKERRQ(ierr);
339*c4762a1bSJed Brown   ierr = TSGetSolveTime(ts, &ftime);CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = TSGetStepNumber(ts,&its);CHKERRQ(ierr);
341*c4762a1bSJed Brown   ierr = PetscPrintf(PETSC_COMM_WORLD, "Number of time steps = %D, final time: %4.2e\nResult:\n\n", its, (double)ftime);CHKERRQ(ierr);
342*c4762a1bSJed Brown   ierr = PrintSolution(x, &user);CHKERRQ(ierr);
343*c4762a1bSJed Brown 
344*c4762a1bSJed Brown   /*
345*c4762a1bSJed Brown    * Free the data structures
346*c4762a1bSJed Brown    */
347*c4762a1bSJed Brown   ierr = VecDestroy(&x);CHKERRQ(ierr);
348*c4762a1bSJed Brown   ierr = MatDestroy(&J);CHKERRQ(ierr);
349*c4762a1bSJed Brown   ierr = TSDestroy(&ts);CHKERRQ(ierr);
350*c4762a1bSJed Brown   ierr = PetscFinalize();
351*c4762a1bSJed Brown   return ierr;
352*c4762a1bSJed Brown }
353*c4762a1bSJed Brown 
354*c4762a1bSJed Brown /*TEST
355*c4762a1bSJed Brown     build:
356*c4762a1bSJed Brown       requires: !single !complex c99
357*c4762a1bSJed Brown 
358*c4762a1bSJed Brown     test:
359*c4762a1bSJed Brown       args: -ts_max_steps 8
360*c4762a1bSJed Brown       output_file: output/ex42.out
361*c4762a1bSJed Brown 
362*c4762a1bSJed Brown TEST*/
363