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