xref: /petsc/src/ts/tutorials/ex20.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\
3c4762a1bSJed Brown Input parameters include:\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /* ------------------------------------------------------------------------
6c4762a1bSJed Brown 
7c4762a1bSJed Brown    This program solves the van der Pol DAE ODE equivalent
8c4762a1bSJed Brown        y' = z                 (1)
9c4762a1bSJed Brown        z' = \mu ((1-y^2)z-y)
10c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
11c4762a1bSJed Brown        y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2),
12c4762a1bSJed Brown    and
13c4762a1bSJed Brown        \mu = 10^6 ( y'(0) ~ -0.6666665432100101).
14c4762a1bSJed Brown    This is a nonlinear equation. The well prepared initial condition gives errors that are not dominated by the first few steps of the method when \mu is large.
15c4762a1bSJed Brown 
16c4762a1bSJed Brown    Notes:
17c4762a1bSJed Brown    This code demonstrates the TS solver interface to an ODE -- RHSFunction for explicit form and IFunction for implicit form.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   ------------------------------------------------------------------------- */
20c4762a1bSJed Brown 
21c4762a1bSJed Brown #include <petscts.h>
22c4762a1bSJed Brown 
23c4762a1bSJed Brown typedef struct _n_User *User;
24c4762a1bSJed Brown struct _n_User {
25c4762a1bSJed Brown   PetscReal mu;
26c4762a1bSJed Brown   PetscReal next_output;
27c4762a1bSJed Brown };
28c4762a1bSJed Brown 
29c4762a1bSJed Brown /*
300e3d61c9SBarry Smith    User-defined routines
31c4762a1bSJed Brown */
329371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx) {
33c4762a1bSJed Brown   User               user = (User)ctx;
34c4762a1bSJed Brown   PetscScalar       *f;
35c4762a1bSJed Brown   const PetscScalar *x;
36c4762a1bSJed Brown 
37c4762a1bSJed Brown   PetscFunctionBeginUser;
389566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
399566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
40c4762a1bSJed Brown   f[0] = x[1];
41c4762a1bSJed Brown   f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0];
429566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
439566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
44c4762a1bSJed Brown   PetscFunctionReturn(0);
45c4762a1bSJed Brown }
46c4762a1bSJed Brown 
479371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
48c4762a1bSJed Brown   User               user = (User)ctx;
49c4762a1bSJed Brown   const PetscScalar *x, *xdot;
50c4762a1bSJed Brown   PetscScalar       *f;
51c4762a1bSJed Brown 
52c4762a1bSJed Brown   PetscFunctionBeginUser;
539566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
549566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
559566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
56c4762a1bSJed Brown   f[0] = xdot[0] - x[1];
57c4762a1bSJed Brown   f[1] = xdot[1] - user->mu * ((1.0 - x[0] * x[0]) * x[1] - x[0]);
589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
599566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
609566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
61c4762a1bSJed Brown   PetscFunctionReturn(0);
62c4762a1bSJed Brown }
63c4762a1bSJed Brown 
649371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
65c4762a1bSJed Brown   User               user     = (User)ctx;
66c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
67c4762a1bSJed Brown   const PetscScalar *x;
68c4762a1bSJed Brown   PetscScalar        J[2][2];
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBeginUser;
719566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
729371c9d4SSatish Balay   J[0][0] = a;
739371c9d4SSatish Balay   J[0][1] = -1.0;
749371c9d4SSatish Balay   J[1][0] = user->mu * (2.0 * x[0] * x[1] + 1.0);
759371c9d4SSatish Balay   J[1][1] = a - user->mu * (1.0 - x[0] * x[0]);
769566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
779566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
809566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
81c4762a1bSJed Brown   if (A != B) {
829566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
839566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
84c4762a1bSJed Brown   }
85c4762a1bSJed Brown   PetscFunctionReturn(0);
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
88c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
899371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) {
90c4762a1bSJed Brown   const PetscScalar *x;
91c4762a1bSJed Brown   PetscReal          tfinal, dt;
92c4762a1bSJed Brown   User               user = (User)ctx;
93c4762a1bSJed Brown   Vec                interpolatedX;
94c4762a1bSJed Brown 
95c4762a1bSJed Brown   PetscFunctionBeginUser;
969566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
979566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
98c4762a1bSJed Brown 
99c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
1009566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &interpolatedX));
1019566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
1029566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX, &x));
1039371c9d4SSatish Balay     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "[%.1f] %" PetscInt_FMT " TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n", (double)user->next_output, step, (double)t, (double)dt, (double)PetscRealPart(x[0]), (double)PetscRealPart(x[1])));
1049566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX, &x));
1059566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
106c4762a1bSJed Brown     user->next_output += 0.1;
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown   PetscFunctionReturn(0);
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
1119371c9d4SSatish Balay int main(int argc, char **argv) {
112c4762a1bSJed Brown   TS             ts; /* nonlinear solver */
113c4762a1bSJed Brown   Vec            x;  /* solution, residual vectors */
114c4762a1bSJed Brown   Mat            A;  /* Jacobian matrix */
115c4762a1bSJed Brown   PetscInt       steps;
116c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
117c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE, implicitform = PETSC_TRUE;
118c4762a1bSJed Brown   PetscScalar   *x_ptr;
119c4762a1bSJed Brown   PetscMPIInt    size;
120c4762a1bSJed Brown   struct _n_User user;
121c4762a1bSJed Brown 
122c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
123c4762a1bSJed Brown      Initialize program
124c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
125327415f7SBarry Smith   PetscFunctionBeginUser;
1269566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1279566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1283c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
131c4762a1bSJed Brown     Set runtime options
132c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
133c4762a1bSJed Brown   user.next_output = 0.0;
134c4762a1bSJed Brown   user.mu          = 1.0e3;
1359566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
1369566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-implicitform", &implicitform, NULL));
137d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Physical parameters", NULL);
1389566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-mu", "Stiffness parameter", "<1.0e6>", user.mu, &user.mu, NULL));
139d0609cedSBarry Smith   PetscOptionsEnd();
140c4762a1bSJed Brown 
141c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
142c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
143c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1449566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1459566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
1469566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1479566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
148c4762a1bSJed Brown 
1499566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
150c4762a1bSJed Brown 
151c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
152c4762a1bSJed Brown      Create timestepping solver context
153c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1549566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
155c4762a1bSJed Brown   if (implicitform) {
1569566063dSJacob Faibussowitsch     PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
1579566063dSJacob Faibussowitsch     PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user));
1589566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSBEULER));
159c4762a1bSJed Brown   } else {
1609566063dSJacob Faibussowitsch     PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
1619566063dSJacob Faibussowitsch     PetscCall(TSSetType(ts, TSRK));
162c4762a1bSJed Brown   }
1639566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
1649566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.001));
1659566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
166*48a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
169c4762a1bSJed Brown      Set initial conditions
170c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1719566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &x_ptr));
172c4762a1bSJed Brown   x_ptr[0] = 2.0;
173c4762a1bSJed Brown   x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
1749566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &x_ptr));
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177c4762a1bSJed Brown      Set runtime options
178c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1799566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
180c4762a1bSJed Brown 
181c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
182c4762a1bSJed Brown      Solve nonlinear system
183c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1849566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
1859566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
1869566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
18763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %" PetscInt_FMT ", ftime %g\n", steps, (double)ftime));
1889566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
189c4762a1bSJed Brown 
190c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
192c4762a1bSJed Brown      are no longer needed.
193c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1949566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
1959566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
1969566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
197c4762a1bSJed Brown 
1989566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
199d0609cedSBarry Smith   return (0);
200c4762a1bSJed Brown }
201c4762a1bSJed Brown 
202c4762a1bSJed Brown /*TEST
203c4762a1bSJed Brown 
204c4762a1bSJed Brown     test:
205c4762a1bSJed Brown       requires: !single
206c4762a1bSJed Brown       args: -mu 1e6
207c4762a1bSJed Brown 
20810b8587eSHendrik Ranocha     test:
20910b8587eSHendrik Ranocha       requires: !single
21010b8587eSHendrik Ranocha       suffix: 2
21110b8587eSHendrik Ranocha       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp
21210b8587eSHendrik Ranocha 
21310b8587eSHendrik Ranocha     test:
21410b8587eSHendrik Ranocha       requires: !single
21510b8587eSHendrik Ranocha       suffix: 3
21610b8587eSHendrik Ranocha       args: -implicitform false -ts_type rk -ts_rk_type 5dp -ts_adapt_type dsp -ts_adapt_dsp_filter H0312
21710b8587eSHendrik Ranocha 
218c4762a1bSJed Brown TEST*/
219