xref: /petsc/src/ts/tutorials/ex42.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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   PetscScalar       *dxdt;
35c4762a1bSJed Brown   const PetscScalar *x;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   PetscFunctionBegin;
38c4762a1bSJed Brown   nb_cells = user->nb_cells;
39c4762a1bSJed Brown   alpha    = user->alpha;
40c4762a1bSJed Brown   beta     = user->beta;
41c4762a1bSJed Brown   rho_a    = user->rho_a;
42c4762a1bSJed Brown   mu_a     = user->mu_a;
43c4762a1bSJed Brown   D_a      = user->D_a;
44c4762a1bSJed Brown   rho_h    = user->rho_h;
45c4762a1bSJed Brown   mu_h     = user->mu_h;
46c4762a1bSJed Brown   D_h      = user->D_h;
47c4762a1bSJed Brown 
489566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
499566063dSJacob Faibussowitsch   PetscCall(VecGetArray(DXDT, &dxdt));
50c4762a1bSJed Brown 
51c4762a1bSJed Brown   for (i = 0 ; i < nb_cells ; i++) {
52c4762a1bSJed Brown     a = x[2*i];
53c4762a1bSJed Brown     h = x[2*i+1];
54c4762a1bSJed Brown     // Reaction:
55c4762a1bSJed Brown     da = alpha * a*a / (1. + beta * h) + rho_a - mu_a * a;
56c4762a1bSJed Brown     dh = alpha * a*a + rho_h - mu_h*h;
57c4762a1bSJed Brown     // Diffusion:
58c4762a1bSJed Brown     d2a = d2h = 0.;
59c4762a1bSJed Brown     if (i > 0) {
60c4762a1bSJed Brown       d2a += (x[2*(i-1)] - a);
61c4762a1bSJed Brown       d2h += (x[2*(i-1)+1] - h);
62c4762a1bSJed Brown     }
63c4762a1bSJed Brown     if (i < nb_cells-1) {
64c4762a1bSJed Brown       d2a += (x[2*(i+1)] - a);
65c4762a1bSJed Brown       d2h += (x[2*(i+1)+1] - h);
66c4762a1bSJed Brown     }
67c4762a1bSJed Brown     dxdt[2*i] = da + D_a*d2a;
68c4762a1bSJed Brown     dxdt[2*i+1] = dh + D_h*d2h;
69c4762a1bSJed Brown   }
709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(DXDT, &dxdt));
719566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
72c4762a1bSJed Brown   PetscFunctionReturn(0);
73c4762a1bSJed Brown }
74c4762a1bSJed Brown 
75c4762a1bSJed Brown PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat J, Mat B, void *ptr)
76c4762a1bSJed Brown {
77c4762a1bSJed Brown   AppCtx            *user = (AppCtx*)ptr;
78c4762a1bSJed Brown   PetscInt          nb_cells, i, idx;
79c4762a1bSJed Brown   PetscReal         alpha, beta;
80c4762a1bSJed Brown   PetscReal         mu_a, D_a;
81c4762a1bSJed Brown   PetscReal         mu_h, D_h;
82c4762a1bSJed Brown   PetscReal         a, h;
83c4762a1bSJed Brown   const PetscScalar *x;
84c4762a1bSJed Brown   PetscScalar       va[4], vh[4];
85c4762a1bSJed Brown   PetscInt          ca[4], ch[4], rowa, rowh;
86c4762a1bSJed Brown 
87c4762a1bSJed Brown   PetscFunctionBegin;
88c4762a1bSJed Brown   nb_cells = user->nb_cells;
89c4762a1bSJed Brown   alpha    = user->alpha;
90c4762a1bSJed Brown   beta     = user->beta;
91c4762a1bSJed Brown   mu_a     = user->mu_a;
92c4762a1bSJed Brown   D_a      = user->D_a;
93c4762a1bSJed Brown   mu_h     = user->mu_h;
94c4762a1bSJed Brown   D_h      = user->D_h;
95c4762a1bSJed Brown 
969566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
97c4762a1bSJed Brown   for (i = 0; i < nb_cells ; ++i) {
98c4762a1bSJed Brown     rowa = 2*i;
99c4762a1bSJed Brown     rowh = 2*i+1;
100c4762a1bSJed Brown     a = x[2*i];
101c4762a1bSJed Brown     h = x[2*i+1];
102c4762a1bSJed Brown     ca[0] = ch[1] = 2*i;
103c4762a1bSJed Brown     va[0] = 2*alpha*a / (1.+beta*h) - mu_a;
104c4762a1bSJed Brown     vh[1] = 2*alpha*a;
105c4762a1bSJed Brown     ca[1] = ch[0] = 2*i+1;
106c4762a1bSJed Brown     va[1] = -beta*alpha*a*a / ((1.+beta*h)*(1.+beta*h));
107c4762a1bSJed Brown     vh[0] = -mu_h;
108c4762a1bSJed Brown     idx = 2;
109c4762a1bSJed Brown     if (i > 0) {
110c4762a1bSJed Brown       ca[idx] = 2*(i-1);
111c4762a1bSJed Brown       ch[idx] = 2*(i-1)+1;
112c4762a1bSJed Brown       va[idx] = D_a;
113c4762a1bSJed Brown       vh[idx] = D_h;
114c4762a1bSJed Brown       va[0] -= D_a;
115c4762a1bSJed Brown       vh[0] -= D_h;
116c4762a1bSJed Brown       idx++;
117c4762a1bSJed Brown     }
118c4762a1bSJed Brown     if (i < nb_cells-1) {
119c4762a1bSJed Brown       ca[idx] = 2*(i+1);
120c4762a1bSJed Brown       ch[idx] = 2*(i+1)+1;
121c4762a1bSJed Brown       va[idx] = D_a;
122c4762a1bSJed Brown       vh[idx] = D_h;
123c4762a1bSJed Brown       va[0] -= D_a;
124c4762a1bSJed Brown       vh[0] -= D_h;
125c4762a1bSJed Brown       idx++;
126c4762a1bSJed Brown     }
1279566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &rowa, idx, ca, va, INSERT_VALUES));
1289566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &rowh, idx, ch, vh, INSERT_VALUES));
129c4762a1bSJed Brown   }
1309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1319566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
1329566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
133c4762a1bSJed Brown   if (J != B) {
1349566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
1359566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
136c4762a1bSJed Brown   }
137c4762a1bSJed Brown   PetscFunctionReturn(0);
138c4762a1bSJed Brown }
139c4762a1bSJed Brown 
140c4762a1bSJed Brown PetscErrorCode DomainErrorFunction(TS ts, PetscReal t, Vec Y, PetscBool *accept)
141c4762a1bSJed Brown {
142c4762a1bSJed Brown   AppCtx            *user;
143c4762a1bSJed Brown   PetscReal         dt;
144c4762a1bSJed Brown   const PetscScalar *x;
145c4762a1bSJed Brown   PetscInt          nb_cells, i;
146c4762a1bSJed Brown 
147362febeeSStefano Zampini   PetscFunctionBegin;
1489566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &user));
149c4762a1bSJed Brown   nb_cells = user->nb_cells;
1509566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &x));
151c4762a1bSJed Brown   for (i = 0 ; i < 2*nb_cells ; ++i) {
152c4762a1bSJed Brown     if (PetscRealPart(x[i]) < 0) {
1539566063dSJacob Faibussowitsch       PetscCall(TSGetTimeStep(ts, &dt));
1549566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ** Domain Error at time %g\n", (double)t));
155c4762a1bSJed Brown       *accept = PETSC_FALSE;
156c4762a1bSJed Brown       break;
157c4762a1bSJed Brown     }
158c4762a1bSJed Brown   }
1599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &x));
160c4762a1bSJed Brown   PetscFunctionReturn(0);
161c4762a1bSJed Brown }
162c4762a1bSJed Brown 
163c4762a1bSJed Brown PetscErrorCode FormInitialState(Vec X, AppCtx* user)
164c4762a1bSJed Brown {
165c4762a1bSJed Brown   PetscRandom    R;
166c4762a1bSJed Brown 
167c4762a1bSJed Brown   PetscFunctionBegin;
1689566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &R));
1699566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(R));
1709566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(R, 0., 10.));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   /*
173c4762a1bSJed Brown    * Initialize the state vector
174c4762a1bSJed Brown    */
1759566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(X, R));
1769566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&R));
177c4762a1bSJed Brown   PetscFunctionReturn(0);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown 
180c4762a1bSJed Brown PetscErrorCode PrintSolution(Vec X, AppCtx *user)
181c4762a1bSJed Brown {
182c4762a1bSJed Brown   const PetscScalar *x;
183c4762a1bSJed Brown   PetscInt          i;
184c4762a1bSJed Brown   PetscInt          nb_cells = user->nb_cells;
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   PetscFunctionBegin;
1879566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1889566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activator,Inhibitor\n"));
189c4762a1bSJed Brown   for (i = 0 ; i < nb_cells ; i++) {
1909566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%5.6e,%5.6e\n", (double)x[2*i], (double)x[2*i+1]));
191c4762a1bSJed Brown   }
1929566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
193c4762a1bSJed Brown   PetscFunctionReturn(0);
194c4762a1bSJed Brown }
195c4762a1bSJed Brown 
196c4762a1bSJed Brown int main(int argc, char **argv)
197c4762a1bSJed Brown {
198c4762a1bSJed Brown   TS             ts;       /* time-stepping context */
199c4762a1bSJed Brown   Vec            x;       /* State vector */
200c4762a1bSJed Brown   Mat            J; /* Jacobian matrix */
201c4762a1bSJed Brown   AppCtx         user; /* user-defined context */
202c4762a1bSJed Brown   PetscReal      ftime;
203c4762a1bSJed Brown   PetscInt       its;
204c4762a1bSJed Brown   PetscMPIInt    size;
205c4762a1bSJed Brown 
206*327415f7SBarry Smith   PetscFunctionBeginUser;
2079566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2089566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2093c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only");
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /*
212c4762a1bSJed Brown    * Allow user to set the grid dimensions and the equations parameters
213c4762a1bSJed Brown    */
214c4762a1bSJed Brown 
215c4762a1bSJed Brown   user.nb_cells = 50;
216c4762a1bSJed Brown   user.alpha = 10.;
217c4762a1bSJed Brown   user.beta = 1.;
218c4762a1bSJed Brown   user.rho_a = 1.;
219c4762a1bSJed Brown   user.rho_h = 2.;
220c4762a1bSJed Brown   user.mu_a = 2.;
221c4762a1bSJed Brown   user.mu_h = 3.;
222c4762a1bSJed Brown   user.D_a = 0.;
223c4762a1bSJed Brown   user.D_h = 30.;
224c4762a1bSJed Brown 
225d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Problem settings", "PROBLEM");
2269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-nb_cells", "Number of cells", "ex42.c",user.nb_cells, &user.nb_cells,NULL));
2279566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-alpha", "Autocatalysis factor", "ex42.c",user.alpha, &user.alpha,NULL));
2289566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-beta", "Inhibition factor", "ex42.c",user.beta, &user.beta,NULL));
2299566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-rho_a", "Default production of the activator", "ex42.c",user.rho_a, &user.rho_a,NULL));
2309566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mu_a", "Degradation rate of the activator", "ex42.c",user.mu_a, &user.mu_a,NULL));
2319566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-D_a", "Diffusion rate of the activator", "ex42.c",user.D_a, &user.D_a,NULL));
2329566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-rho_h", "Default production of the inhibitor", "ex42.c",user.rho_h, &user.rho_h,NULL));
2339566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mu_h", "Degradation rate of the inhibitor", "ex42.c",user.mu_h, &user.mu_h,NULL));
2349566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-D_h", "Diffusion rate of the inhibitor", "ex42.c",user.D_h, &user.D_h,NULL));
235d0609cedSBarry Smith   PetscOptionsEnd();
236c4762a1bSJed Brown 
23763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "nb_cells: %" PetscInt_FMT "\n", user.nb_cells));
2389566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "alpha: %5.5g\n", (double)user.alpha));
2399566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "beta:  %5.5g\n", (double)user.beta));
2409566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "rho_a: %5.5g\n", (double)user.rho_a));
2419566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu_a:  %5.5g\n", (double)user.mu_a));
2429566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "D_a:   %5.5g\n", (double)user.D_a));
2439566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "rho_h: %5.5g\n", (double)user.rho_h));
2449566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu_h:  %5.5g\n", (double)user.mu_h));
2459566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "D_h:   %5.5g\n", (double)user.D_h));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown   /*
248c4762a1bSJed Brown    * Create vector to hold the solution
249c4762a1bSJed Brown    */
2509566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 2*user.nb_cells, &x));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /*
253c4762a1bSJed Brown    * Create time-stepper context
254c4762a1bSJed Brown    */
2559566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2569566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
257c4762a1bSJed Brown 
258c4762a1bSJed Brown   /*
259c4762a1bSJed Brown    * Tell the time-stepper context where to compute the solution
260c4762a1bSJed Brown    */
2619566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
262c4762a1bSJed Brown 
263c4762a1bSJed Brown   /*
264c4762a1bSJed Brown    * Allocate the jacobian matrix
265c4762a1bSJed Brown    */
2669566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 2*user.nb_cells, 2*user.nb_cells, 4, 0, &J));
267c4762a1bSJed Brown 
268c4762a1bSJed Brown   /*
269c4762a1bSJed Brown    * Provide the call-back for the non-linear function we are evaluating.
270c4762a1bSJed Brown    */
2719566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
272c4762a1bSJed Brown 
273c4762a1bSJed Brown   /*
274c4762a1bSJed Brown    * Set the Jacobian matrix and the function user to compute Jacobians
275c4762a1bSJed Brown    */
2769566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, &user));
277c4762a1bSJed Brown 
278c4762a1bSJed Brown   /*
279c4762a1bSJed Brown    * Set the function checking the domain
280c4762a1bSJed Brown    */
2819566063dSJacob Faibussowitsch   PetscCall(TSSetFunctionDomainError(ts, &DomainErrorFunction));
282c4762a1bSJed Brown 
283c4762a1bSJed Brown   /*
284c4762a1bSJed Brown    * Initialize the problem with random values
285c4762a1bSJed Brown    */
2869566063dSJacob Faibussowitsch   PetscCall(FormInitialState(x, &user));
287c4762a1bSJed Brown 
288c4762a1bSJed Brown   /*
289c4762a1bSJed Brown    * Read the solver type from options
290c4762a1bSJed Brown    */
2919566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSPSEUDO));
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   /*
294c4762a1bSJed Brown    * Set a large number of timesteps and final duration time to insure
295c4762a1bSJed Brown    * convergenge to steady state
296c4762a1bSJed Brown    */
2979566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 2147483647));
2989566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1.e12));
2999566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER));
300c4762a1bSJed Brown 
301c4762a1bSJed Brown   /*
302c4762a1bSJed Brown    * Set a larger number of potential errors
303c4762a1bSJed Brown    */
3049566063dSJacob Faibussowitsch   PetscCall(TSSetMaxStepRejections(ts, 50));
305c4762a1bSJed Brown 
306c4762a1bSJed Brown   /*
307c4762a1bSJed Brown    * Also start with a very small dt
308c4762a1bSJed Brown    */
3099566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.05));
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   /*
312c4762a1bSJed Brown    * Set a larger time step increment
313c4762a1bSJed Brown    */
3149566063dSJacob Faibussowitsch   PetscCall(TSPseudoSetTimeStepIncrement(ts, 1.5));
315c4762a1bSJed Brown 
316c4762a1bSJed Brown   /*
317c4762a1bSJed Brown    * Let the user personalise TS
318c4762a1bSJed Brown    */
3199566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
320c4762a1bSJed Brown 
321c4762a1bSJed Brown   /*
322c4762a1bSJed Brown    * Set the context for the time stepper
323c4762a1bSJed Brown    */
3249566063dSJacob Faibussowitsch   PetscCall(TSSetApplicationContext(ts, &user));
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   /*
327c4762a1bSJed Brown    * Setup the time stepper, ready for evaluation
328c4762a1bSJed Brown    */
3299566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
330c4762a1bSJed Brown 
331c4762a1bSJed Brown   /*
332c4762a1bSJed Brown    * Perform the solve.
333c4762a1bSJed Brown    */
3349566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
3359566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
3369566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&its));
33763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of time steps = %" PetscInt_FMT ", final time: %4.2e\nResult:\n\n", its, (double)ftime));
3389566063dSJacob Faibussowitsch   PetscCall(PrintSolution(x, &user));
339c4762a1bSJed Brown 
340c4762a1bSJed Brown   /*
341c4762a1bSJed Brown    * Free the data structures
342c4762a1bSJed Brown    */
3439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3449566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
3459566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3469566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
347b122ec5aSJacob Faibussowitsch   return 0;
348c4762a1bSJed Brown }
349c4762a1bSJed Brown 
350c4762a1bSJed Brown /*TEST
351c4762a1bSJed Brown     build:
352f56ea12dSJed Brown       requires: !single !complex
353c4762a1bSJed Brown 
354c4762a1bSJed Brown     test:
355c4762a1bSJed Brown       args: -ts_max_steps 8
356c4762a1bSJed Brown       output_file: output/ex42.out
357c4762a1bSJed Brown 
358c4762a1bSJed Brown TEST*/
359