xref: /petsc/src/ts/tutorials/ex42.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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 
269371c9d4SSatish Balay PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec DXDT, void *ptr) {
27c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
28c4762a1bSJed Brown   PetscInt           nb_cells, i;
29c4762a1bSJed Brown   PetscReal          alpha, beta;
30c4762a1bSJed Brown   PetscReal          rho_a, mu_a, D_a;
31c4762a1bSJed Brown   PetscReal          rho_h, mu_h, D_h;
32c4762a1bSJed Brown   PetscReal          a, h, da, dh, d2a, d2h;
33c4762a1bSJed Brown   PetscScalar       *dxdt;
34c4762a1bSJed Brown   const PetscScalar *x;
35c4762a1bSJed Brown 
367510d9b0SBarry Smith   PetscFunctionBeginUser;
37c4762a1bSJed Brown   nb_cells = user->nb_cells;
38c4762a1bSJed Brown   alpha    = user->alpha;
39c4762a1bSJed Brown   beta     = user->beta;
40c4762a1bSJed Brown   rho_a    = user->rho_a;
41c4762a1bSJed Brown   mu_a     = user->mu_a;
42c4762a1bSJed Brown   D_a      = user->D_a;
43c4762a1bSJed Brown   rho_h    = user->rho_h;
44c4762a1bSJed Brown   mu_h     = user->mu_h;
45c4762a1bSJed Brown   D_h      = user->D_h;
46c4762a1bSJed Brown 
479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
489566063dSJacob Faibussowitsch   PetscCall(VecGetArray(DXDT, &dxdt));
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   for (i = 0; i < nb_cells; i++) {
51c4762a1bSJed Brown     a   = x[2 * i];
52c4762a1bSJed Brown     h   = x[2 * i + 1];
53c4762a1bSJed Brown     // Reaction:
54c4762a1bSJed Brown     da  = alpha * a * a / (1. + beta * h) + rho_a - mu_a * a;
55c4762a1bSJed Brown     dh  = alpha * a * a + rho_h - mu_h * h;
56c4762a1bSJed Brown     // Diffusion:
57c4762a1bSJed Brown     d2a = d2h = 0.;
58c4762a1bSJed Brown     if (i > 0) {
59c4762a1bSJed Brown       d2a += (x[2 * (i - 1)] - a);
60c4762a1bSJed Brown       d2h += (x[2 * (i - 1) + 1] - h);
61c4762a1bSJed Brown     }
62c4762a1bSJed Brown     if (i < nb_cells - 1) {
63c4762a1bSJed Brown       d2a += (x[2 * (i + 1)] - a);
64c4762a1bSJed Brown       d2h += (x[2 * (i + 1) + 1] - h);
65c4762a1bSJed Brown     }
66c4762a1bSJed Brown     dxdt[2 * i]     = da + D_a * d2a;
67c4762a1bSJed Brown     dxdt[2 * i + 1] = dh + D_h * d2h;
68c4762a1bSJed Brown   }
699566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(DXDT, &dxdt));
709566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
71c4762a1bSJed Brown   PetscFunctionReturn(0);
72c4762a1bSJed Brown }
73c4762a1bSJed Brown 
749371c9d4SSatish Balay PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat J, Mat B, void *ptr) {
75c4762a1bSJed Brown   AppCtx            *user = (AppCtx *)ptr;
76c4762a1bSJed Brown   PetscInt           nb_cells, i, idx;
77c4762a1bSJed Brown   PetscReal          alpha, beta;
78c4762a1bSJed Brown   PetscReal          mu_a, D_a;
79c4762a1bSJed Brown   PetscReal          mu_h, D_h;
80c4762a1bSJed Brown   PetscReal          a, h;
81c4762a1bSJed Brown   const PetscScalar *x;
82c4762a1bSJed Brown   PetscScalar        va[4], vh[4];
83c4762a1bSJed Brown   PetscInt           ca[4], ch[4], rowa, rowh;
84c4762a1bSJed Brown 
857510d9b0SBarry Smith   PetscFunctionBeginUser;
86c4762a1bSJed Brown   nb_cells = user->nb_cells;
87c4762a1bSJed Brown   alpha    = user->alpha;
88c4762a1bSJed Brown   beta     = user->beta;
89c4762a1bSJed Brown   mu_a     = user->mu_a;
90c4762a1bSJed Brown   D_a      = user->D_a;
91c4762a1bSJed Brown   mu_h     = user->mu_h;
92c4762a1bSJed Brown   D_h      = user->D_h;
93c4762a1bSJed Brown 
949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
95c4762a1bSJed Brown   for (i = 0; i < nb_cells; ++i) {
96c4762a1bSJed Brown     rowa  = 2 * i;
97c4762a1bSJed Brown     rowh  = 2 * i + 1;
98c4762a1bSJed Brown     a     = x[2 * i];
99c4762a1bSJed Brown     h     = x[2 * i + 1];
100c4762a1bSJed Brown     ca[0] = ch[1] = 2 * i;
101c4762a1bSJed Brown     va[0]         = 2 * alpha * a / (1. + beta * h) - mu_a;
102c4762a1bSJed Brown     vh[1]         = 2 * alpha * a;
103c4762a1bSJed Brown     ca[1] = ch[0] = 2 * i + 1;
104c4762a1bSJed Brown     va[1]         = -beta * alpha * a * a / ((1. + beta * h) * (1. + beta * h));
105c4762a1bSJed Brown     vh[0]         = -mu_h;
106c4762a1bSJed Brown     idx           = 2;
107c4762a1bSJed Brown     if (i > 0) {
108c4762a1bSJed Brown       ca[idx] = 2 * (i - 1);
109c4762a1bSJed Brown       ch[idx] = 2 * (i - 1) + 1;
110c4762a1bSJed Brown       va[idx] = D_a;
111c4762a1bSJed Brown       vh[idx] = D_h;
112c4762a1bSJed Brown       va[0] -= D_a;
113c4762a1bSJed Brown       vh[0] -= D_h;
114c4762a1bSJed Brown       idx++;
115c4762a1bSJed Brown     }
116c4762a1bSJed Brown     if (i < nb_cells - 1) {
117c4762a1bSJed Brown       ca[idx] = 2 * (i + 1);
118c4762a1bSJed Brown       ch[idx] = 2 * (i + 1) + 1;
119c4762a1bSJed Brown       va[idx] = D_a;
120c4762a1bSJed Brown       vh[idx] = D_h;
121c4762a1bSJed Brown       va[0] -= D_a;
122c4762a1bSJed Brown       vh[0] -= D_h;
123c4762a1bSJed Brown       idx++;
124c4762a1bSJed Brown     }
1259566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &rowa, idx, ca, va, INSERT_VALUES));
1269566063dSJacob Faibussowitsch     PetscCall(MatSetValues(B, 1, &rowh, idx, ch, vh, INSERT_VALUES));
127c4762a1bSJed Brown   }
1289566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1309566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
131c4762a1bSJed Brown   if (J != B) {
1329566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
1339566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
134c4762a1bSJed Brown   }
135c4762a1bSJed Brown   PetscFunctionReturn(0);
136c4762a1bSJed Brown }
137c4762a1bSJed Brown 
1389371c9d4SSatish Balay PetscErrorCode DomainErrorFunction(TS ts, PetscReal t, Vec Y, PetscBool *accept) {
139c4762a1bSJed Brown   AppCtx            *user;
140c4762a1bSJed Brown   PetscReal          dt;
141c4762a1bSJed Brown   const PetscScalar *x;
142c4762a1bSJed Brown   PetscInt           nb_cells, i;
143c4762a1bSJed Brown 
1447510d9b0SBarry Smith   PetscFunctionBeginUser;
1459566063dSJacob Faibussowitsch   PetscCall(TSGetApplicationContext(ts, &user));
146c4762a1bSJed Brown   nb_cells = user->nb_cells;
1479566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Y, &x));
148c4762a1bSJed Brown   for (i = 0; i < 2 * nb_cells; ++i) {
149c4762a1bSJed Brown     if (PetscRealPart(x[i]) < 0) {
1509566063dSJacob Faibussowitsch       PetscCall(TSGetTimeStep(ts, &dt));
1519566063dSJacob Faibussowitsch       PetscCall(PetscPrintf(PETSC_COMM_WORLD, " ** Domain Error at time %g\n", (double)t));
152c4762a1bSJed Brown       *accept = PETSC_FALSE;
153c4762a1bSJed Brown       break;
154c4762a1bSJed Brown     }
155c4762a1bSJed Brown   }
1569566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Y, &x));
157c4762a1bSJed Brown   PetscFunctionReturn(0);
158c4762a1bSJed Brown }
159c4762a1bSJed Brown 
1609371c9d4SSatish Balay PetscErrorCode FormInitialState(Vec X, AppCtx *user) {
161c4762a1bSJed Brown   PetscRandom R;
162c4762a1bSJed Brown 
1637510d9b0SBarry Smith   PetscFunctionBeginUser;
1649566063dSJacob Faibussowitsch   PetscCall(PetscRandomCreate(PETSC_COMM_WORLD, &R));
1659566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetFromOptions(R));
1669566063dSJacob Faibussowitsch   PetscCall(PetscRandomSetInterval(R, 0., 10.));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /*
169c4762a1bSJed Brown    * Initialize the state vector
170c4762a1bSJed Brown    */
1719566063dSJacob Faibussowitsch   PetscCall(VecSetRandom(X, R));
1729566063dSJacob Faibussowitsch   PetscCall(PetscRandomDestroy(&R));
173c4762a1bSJed Brown   PetscFunctionReturn(0);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
1769371c9d4SSatish Balay PetscErrorCode PrintSolution(Vec X, AppCtx *user) {
177c4762a1bSJed Brown   const PetscScalar *x;
178c4762a1bSJed Brown   PetscInt           i;
179c4762a1bSJed Brown   PetscInt           nb_cells = user->nb_cells;
180c4762a1bSJed Brown 
1817510d9b0SBarry Smith   PetscFunctionBeginUser;
1829566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1839566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Activator,Inhibitor\n"));
184*48a46eb9SPierre Jolivet   for (i = 0; i < nb_cells; i++) PetscCall(PetscPrintf(PETSC_COMM_WORLD, "%5.6e,%5.6e\n", (double)x[2 * i], (double)x[2 * i + 1]));
1859566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
186c4762a1bSJed Brown   PetscFunctionReturn(0);
187c4762a1bSJed Brown }
188c4762a1bSJed Brown 
1899371c9d4SSatish Balay int main(int argc, char **argv) {
190c4762a1bSJed Brown   TS          ts;   /* time-stepping context */
191c4762a1bSJed Brown   Vec         x;    /* State vector */
192c4762a1bSJed Brown   Mat         J;    /* Jacobian matrix */
193c4762a1bSJed Brown   AppCtx      user; /* user-defined context */
194c4762a1bSJed Brown   PetscReal   ftime;
195c4762a1bSJed Brown   PetscInt    its;
196c4762a1bSJed Brown   PetscMPIInt size;
197c4762a1bSJed Brown 
198327415f7SBarry Smith   PetscFunctionBeginUser;
1999566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
2009566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
2013c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only");
202c4762a1bSJed Brown 
203c4762a1bSJed Brown   /*
204c4762a1bSJed Brown    * Allow user to set the grid dimensions and the equations parameters
205c4762a1bSJed Brown    */
206c4762a1bSJed Brown 
207c4762a1bSJed Brown   user.nb_cells = 50;
208c4762a1bSJed Brown   user.alpha    = 10.;
209c4762a1bSJed Brown   user.beta     = 1.;
210c4762a1bSJed Brown   user.rho_a    = 1.;
211c4762a1bSJed Brown   user.rho_h    = 2.;
212c4762a1bSJed Brown   user.mu_a     = 2.;
213c4762a1bSJed Brown   user.mu_h     = 3.;
214c4762a1bSJed Brown   user.D_a      = 0.;
215c4762a1bSJed Brown   user.D_h      = 30.;
216c4762a1bSJed Brown 
217d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, "", "Problem settings", "PROBLEM");
2189566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-nb_cells", "Number of cells", "ex42.c", user.nb_cells, &user.nb_cells, NULL));
2199566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-alpha", "Autocatalysis factor", "ex42.c", user.alpha, &user.alpha, NULL));
2209566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-beta", "Inhibition factor", "ex42.c", user.beta, &user.beta, NULL));
2219566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-rho_a", "Default production of the activator", "ex42.c", user.rho_a, &user.rho_a, NULL));
2229566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mu_a", "Degradation rate of the activator", "ex42.c", user.mu_a, &user.mu_a, NULL));
2239566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-D_a", "Diffusion rate of the activator", "ex42.c", user.D_a, &user.D_a, NULL));
2249566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-rho_h", "Default production of the inhibitor", "ex42.c", user.rho_h, &user.rho_h, NULL));
2259566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mu_h", "Degradation rate of the inhibitor", "ex42.c", user.mu_h, &user.mu_h, NULL));
2269566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-D_h", "Diffusion rate of the inhibitor", "ex42.c", user.D_h, &user.D_h, NULL));
227d0609cedSBarry Smith   PetscOptionsEnd();
228c4762a1bSJed Brown 
22963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "nb_cells: %" PetscInt_FMT "\n", user.nb_cells));
2309566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "alpha: %5.5g\n", (double)user.alpha));
2319566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "beta:  %5.5g\n", (double)user.beta));
2329566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "rho_a: %5.5g\n", (double)user.rho_a));
2339566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu_a:  %5.5g\n", (double)user.mu_a));
2349566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "D_a:   %5.5g\n", (double)user.D_a));
2359566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "rho_h: %5.5g\n", (double)user.rho_h));
2369566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu_h:  %5.5g\n", (double)user.mu_h));
2379566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "D_h:   %5.5g\n", (double)user.D_h));
238c4762a1bSJed Brown 
239c4762a1bSJed Brown   /*
240c4762a1bSJed Brown    * Create vector to hold the solution
241c4762a1bSJed Brown    */
2429566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 2 * user.nb_cells, &x));
243c4762a1bSJed Brown 
244c4762a1bSJed Brown   /*
245c4762a1bSJed Brown    * Create time-stepper context
246c4762a1bSJed Brown    */
2479566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2489566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
249c4762a1bSJed Brown 
250c4762a1bSJed Brown   /*
251c4762a1bSJed Brown    * Tell the time-stepper context where to compute the solution
252c4762a1bSJed Brown    */
2539566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, x));
254c4762a1bSJed Brown 
255c4762a1bSJed Brown   /*
256c4762a1bSJed Brown    * Allocate the jacobian matrix
257c4762a1bSJed Brown    */
2589566063dSJacob Faibussowitsch   PetscCall(MatCreateSeqAIJ(PETSC_COMM_WORLD, 2 * user.nb_cells, 2 * user.nb_cells, 4, 0, &J));
259c4762a1bSJed Brown 
260c4762a1bSJed Brown   /*
261c4762a1bSJed Brown    * Provide the call-back for the non-linear function we are evaluating.
262c4762a1bSJed Brown    */
2639566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /*
266c4762a1bSJed Brown    * Set the Jacobian matrix and the function user to compute Jacobians
267c4762a1bSJed Brown    */
2689566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobian(ts, J, J, RHSJacobian, &user));
269c4762a1bSJed Brown 
270c4762a1bSJed Brown   /*
271c4762a1bSJed Brown    * Set the function checking the domain
272c4762a1bSJed Brown    */
2739566063dSJacob Faibussowitsch   PetscCall(TSSetFunctionDomainError(ts, &DomainErrorFunction));
274c4762a1bSJed Brown 
275c4762a1bSJed Brown   /*
276c4762a1bSJed Brown    * Initialize the problem with random values
277c4762a1bSJed Brown    */
2789566063dSJacob Faibussowitsch   PetscCall(FormInitialState(x, &user));
279c4762a1bSJed Brown 
280c4762a1bSJed Brown   /*
281c4762a1bSJed Brown    * Read the solver type from options
282c4762a1bSJed Brown    */
2839566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSPSEUDO));
284c4762a1bSJed Brown 
285c4762a1bSJed Brown   /*
286c4762a1bSJed Brown    * Set a large number of timesteps and final duration time to insure
287c4762a1bSJed Brown    * convergenge to steady state
288c4762a1bSJed Brown    */
2899566063dSJacob Faibussowitsch   PetscCall(TSSetMaxSteps(ts, 2147483647));
2909566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1.e12));
2919566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
292c4762a1bSJed Brown 
293c4762a1bSJed Brown   /*
294c4762a1bSJed Brown    * Set a larger number of potential errors
295c4762a1bSJed Brown    */
2969566063dSJacob Faibussowitsch   PetscCall(TSSetMaxStepRejections(ts, 50));
297c4762a1bSJed Brown 
298c4762a1bSJed Brown   /*
299c4762a1bSJed Brown    * Also start with a very small dt
300c4762a1bSJed Brown    */
3019566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.05));
302c4762a1bSJed Brown 
303c4762a1bSJed Brown   /*
304c4762a1bSJed Brown    * Set a larger time step increment
305c4762a1bSJed Brown    */
3069566063dSJacob Faibussowitsch   PetscCall(TSPseudoSetTimeStepIncrement(ts, 1.5));
307c4762a1bSJed Brown 
308c4762a1bSJed Brown   /*
309c4762a1bSJed Brown    * Let the user personalise TS
310c4762a1bSJed Brown    */
3119566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
312c4762a1bSJed Brown 
313c4762a1bSJed Brown   /*
314c4762a1bSJed Brown    * Set the context for the time stepper
315c4762a1bSJed Brown    */
3169566063dSJacob Faibussowitsch   PetscCall(TSSetApplicationContext(ts, &user));
317c4762a1bSJed Brown 
318c4762a1bSJed Brown   /*
319c4762a1bSJed Brown    * Setup the time stepper, ready for evaluation
320c4762a1bSJed Brown    */
3219566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
322c4762a1bSJed Brown 
323c4762a1bSJed Brown   /*
324c4762a1bSJed Brown    * Perform the solve.
325c4762a1bSJed Brown    */
3269566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
3279566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
3289566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &its));
32963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Number of time steps = %" PetscInt_FMT ", final time: %4.2e\nResult:\n\n", its, (double)ftime));
3309566063dSJacob Faibussowitsch   PetscCall(PrintSolution(x, &user));
331c4762a1bSJed Brown 
332c4762a1bSJed Brown   /*
333c4762a1bSJed Brown    * Free the data structures
334c4762a1bSJed Brown    */
3359566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
3369566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&J));
3379566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
3389566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
339b122ec5aSJacob Faibussowitsch   return 0;
340c4762a1bSJed Brown }
341c4762a1bSJed Brown 
342c4762a1bSJed Brown /*TEST
343c4762a1bSJed Brown     build:
344f56ea12dSJed Brown       requires: !single !complex
345c4762a1bSJed Brown 
346c4762a1bSJed Brown     test:
347c4762a1bSJed Brown       args: -ts_max_steps 8
348c4762a1bSJed Brown       output_file: output/ex42.out
349c4762a1bSJed Brown 
350c4762a1bSJed Brown TEST*/
351