xref: /petsc/src/ts/tutorials/ex16.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation and demonstrate IMEX.\n\
3c4762a1bSJed Brown Input parameters include:\n\
4c4762a1bSJed Brown       -mu : stiffness parameter\n\n";
5c4762a1bSJed Brown 
6c4762a1bSJed Brown /* ------------------------------------------------------------------------
7c4762a1bSJed Brown 
8c4762a1bSJed Brown    This program solves the van der Pol equation
9c4762a1bSJed Brown        y'' - \mu ((1-y^2)*y' - y) = 0        (1)
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    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.
13c4762a1bSJed Brown 
14c4762a1bSJed Brown    Notes:
15c4762a1bSJed Brown    This code demonstrates the TS solver interface to two variants of
16c4762a1bSJed Brown    linear problems, u_t = f(u,t), namely turning (1) into a system of
17c4762a1bSJed Brown    first order differential equations,
18c4762a1bSJed Brown 
19c4762a1bSJed Brown    [ y' ] = [          z            ]
20c4762a1bSJed Brown    [ z' ]   [ \mu ((1 - y^2) z - y) ]
21c4762a1bSJed Brown 
22c4762a1bSJed Brown    which then we can write as a vector equation
23c4762a1bSJed Brown 
24c4762a1bSJed Brown    [ u_1' ] = [             u_2           ]  (2)
25c4762a1bSJed Brown    [ u_2' ]   [ \mu (1 - u_1^2) u_2 - u_1 ]
26c4762a1bSJed Brown 
27c4762a1bSJed Brown    which is now in the desired form of u_t = f(u,t). One way that we
28c4762a1bSJed Brown    can split f(u,t) in (2) is to split by component,
29c4762a1bSJed Brown 
30c4762a1bSJed Brown    [ u_1' ] = [ u_2 ] + [            0                ]
31c4762a1bSJed Brown    [ u_2' ]   [  0  ]   [ \mu ((1 - u_1^2) u_2 - u_1) ]
32c4762a1bSJed Brown 
33c4762a1bSJed Brown    where
34c4762a1bSJed Brown 
355ab1ac2bSHong Zhang    [ G(u,t) ] = [ u_2 ]
36c4762a1bSJed Brown                 [  0  ]
37c4762a1bSJed Brown 
38c4762a1bSJed Brown    and
39c4762a1bSJed Brown 
405ab1ac2bSHong Zhang    [ F(u',u,t) ] = [ u_1' ] - [            0                ]
41c4762a1bSJed Brown                    [ u_2' ]   [ \mu ((1 - u_1^2) u_2 - u_1) ]
42c4762a1bSJed Brown 
435ab1ac2bSHong Zhang    Using the definition of the Jacobian of F (from the PETSc user manual),
445ab1ac2bSHong Zhang    in the equation F(u',u,t) = G(u,t),
45c4762a1bSJed Brown 
465ab1ac2bSHong Zhang               dF   dF
475ab1ac2bSHong Zhang    J(F) = a * -- - --
48c4762a1bSJed Brown               du'  du
49c4762a1bSJed Brown 
50c4762a1bSJed Brown    where d is the partial derivative. In this example,
51c4762a1bSJed Brown 
525ab1ac2bSHong Zhang    dF   [ 1 ; 0 ]
53c4762a1bSJed Brown    -- = [       ]
54c4762a1bSJed Brown    du'  [ 0 ; 1 ]
55c4762a1bSJed Brown 
565ab1ac2bSHong Zhang    dF   [       0             ;         0        ]
57c4762a1bSJed Brown    -- = [                                        ]
58c4762a1bSJed Brown    du   [ -\mu (2*u_1*u_2 + 1);  \mu (1 - u_1^2) ]
59c4762a1bSJed Brown 
60c4762a1bSJed Brown    Hence,
61c4762a1bSJed Brown 
62c4762a1bSJed Brown           [      a             ;          0          ]
635ab1ac2bSHong Zhang    J(F) = [                                          ]
64c4762a1bSJed Brown           [ \mu (2*u_1*u_2 + 1); a - \mu (1 - u_1^2) ]
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   ------------------------------------------------------------------------- */
67c4762a1bSJed Brown 
68c4762a1bSJed Brown #include <petscts.h>
69c4762a1bSJed Brown 
70c4762a1bSJed Brown typedef struct _n_User *User;
71c4762a1bSJed Brown struct _n_User {
72c4762a1bSJed Brown   PetscReal mu;
73c4762a1bSJed Brown   PetscBool imex;
74c4762a1bSJed Brown   PetscReal next_output;
75c4762a1bSJed Brown };
76c4762a1bSJed Brown 
77c4762a1bSJed Brown /*
780e3d61c9SBarry Smith    User-defined routines
79c4762a1bSJed Brown */
809371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx) {
81c4762a1bSJed Brown   User               user = (User)ctx;
82c4762a1bSJed Brown   PetscScalar       *f;
83c4762a1bSJed Brown   const PetscScalar *x;
84c4762a1bSJed Brown 
85c4762a1bSJed Brown   PetscFunctionBeginUser;
869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
879566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
88c4762a1bSJed Brown   f[0] = (user->imex ? x[1] : 0);
89c4762a1bSJed Brown   f[1] = 0.0;
909566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
919566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
92c4762a1bSJed Brown   PetscFunctionReturn(0);
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
959371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
96c4762a1bSJed Brown   User               user = (User)ctx;
97c4762a1bSJed Brown   const PetscScalar *x, *xdot;
98c4762a1bSJed Brown   PetscScalar       *f;
99c4762a1bSJed Brown 
100c4762a1bSJed Brown   PetscFunctionBeginUser;
1019566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1029566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
1039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
104c4762a1bSJed Brown   f[0] = xdot[0] + (user->imex ? 0 : x[1]);
105c4762a1bSJed Brown   f[1] = xdot[1] - user->mu * ((1. - x[0] * x[0]) * x[1] - x[0]);
1069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
1079566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
1089566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
109c4762a1bSJed Brown   PetscFunctionReturn(0);
110c4762a1bSJed Brown }
111c4762a1bSJed Brown 
1129371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
113c4762a1bSJed Brown   User               user     = (User)ctx;
114c4762a1bSJed Brown   PetscReal          mu       = user->mu;
115c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
116c4762a1bSJed Brown   const PetscScalar *x;
117c4762a1bSJed Brown   PetscScalar        J[2][2];
118c4762a1bSJed Brown 
119c4762a1bSJed Brown   PetscFunctionBeginUser;
1209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
1219371c9d4SSatish Balay   J[0][0] = a;
1229371c9d4SSatish Balay   J[0][1] = (user->imex ? 0 : 1.);
1239371c9d4SSatish Balay   J[1][0] = mu * (2. * x[0] * x[1] + 1.);
1249371c9d4SSatish Balay   J[1][1] = a - mu * (1. - x[0] * x[0]);
1259566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
1269566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
127c4762a1bSJed Brown 
1289566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1299566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
130c4762a1bSJed Brown   if (A != B) {
1319566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1329566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
133c4762a1bSJed Brown   }
134c4762a1bSJed Brown   PetscFunctionReturn(0);
135c4762a1bSJed Brown }
136c4762a1bSJed Brown 
1379371c9d4SSatish Balay static PetscErrorCode RegisterMyARK2(void) {
138c4762a1bSJed Brown   PetscFunctionBeginUser;
139c4762a1bSJed Brown   {
1409371c9d4SSatish Balay     const PetscReal A[3][3] =
1419371c9d4SSatish Balay       {
1429371c9d4SSatish Balay         {0,                      0,    0},
143c4762a1bSJed Brown         {0.41421356237309504880, 0,    0},
1449371c9d4SSatish Balay         {0.75,                   0.25, 0}
1459371c9d4SSatish Balay     },
1469371c9d4SSatish Balay                     At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL;
1479566063dSJacob Faibussowitsch     PetscCall(TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL));
148c4762a1bSJed Brown   }
149c4762a1bSJed Brown   PetscFunctionReturn(0);
150c4762a1bSJed Brown }
151c4762a1bSJed Brown 
152c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
1539371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) {
154c4762a1bSJed Brown   const PetscScalar *x;
155c4762a1bSJed Brown   PetscReal          tfinal, dt;
156c4762a1bSJed Brown   User               user = (User)ctx;
157c4762a1bSJed Brown   Vec                interpolatedX;
158c4762a1bSJed Brown 
159c4762a1bSJed Brown   PetscFunctionBeginUser;
1609566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1619566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
1649566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &interpolatedX));
1659566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
1669566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX, &x));
16763a3b9bcSJacob Faibussowitsch     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])));
1689566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX, &x));
1699566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
170c4762a1bSJed Brown 
171c4762a1bSJed Brown     user->next_output += 0.1;
172c4762a1bSJed Brown   }
173c4762a1bSJed Brown   PetscFunctionReturn(0);
174c4762a1bSJed Brown }
175c4762a1bSJed Brown 
1769371c9d4SSatish Balay int main(int argc, char **argv) {
177c4762a1bSJed Brown   TS             ts; /* nonlinear solver */
178c4762a1bSJed Brown   Vec            x;  /* solution, residual vectors */
179c4762a1bSJed Brown   Mat            A;  /* Jacobian matrix */
180c4762a1bSJed Brown   PetscInt       steps;
181c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
182c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
183c4762a1bSJed Brown   PetscScalar   *x_ptr;
184c4762a1bSJed Brown   PetscMPIInt    size;
185c4762a1bSJed Brown   struct _n_User user;
186c4762a1bSJed Brown 
187c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188c4762a1bSJed Brown      Initialize program
189c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
190327415f7SBarry Smith   PetscFunctionBeginUser;
1919566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1929566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1933c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
194c4762a1bSJed Brown 
1959566063dSJacob Faibussowitsch   PetscCall(RegisterMyARK2());
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown     Set runtime options
199c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200c4762a1bSJed Brown   user.mu          = 1000.0;
201c4762a1bSJed Brown   user.imex        = PETSC_TRUE;
202c4762a1bSJed Brown   user.next_output = 0.0;
203c4762a1bSJed Brown 
2049566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
2059566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-imex", &user.imex, NULL));
2069566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
207c4762a1bSJed Brown 
208c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
210c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2119566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2129566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
2139566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
2149566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
2159566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
216c4762a1bSJed Brown 
217c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218c4762a1bSJed Brown      Create timestepping solver context
219c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2209566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2219566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
2229566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
2239566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
2249566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user));
2259566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
2269566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
227*48a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
228c4762a1bSJed Brown 
229c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
230c4762a1bSJed Brown      Set initial conditions
231c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2329566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &x_ptr));
233c4762a1bSJed Brown   x_ptr[0] = 2.0;
234c4762a1bSJed Brown   x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
2359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &x_ptr));
2369566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.01));
237c4762a1bSJed Brown 
238c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239c4762a1bSJed Brown      Set runtime options
240c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2419566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
242c4762a1bSJed Brown 
243c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown      Solve nonlinear system
245c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2469566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
2479566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2489566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
24963a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime));
2509566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
251c4762a1bSJed Brown 
252c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
253c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
254c4762a1bSJed Brown      are no longer needed.
255c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2569566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2579566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2589566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
259c4762a1bSJed Brown 
2609566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
261b122ec5aSJacob Faibussowitsch   return 0;
262c4762a1bSJed Brown }
263c4762a1bSJed Brown 
264c4762a1bSJed Brown /*TEST
265c4762a1bSJed Brown 
266c4762a1bSJed Brown     test:
267c4762a1bSJed Brown       args: -ts_type arkimex -ts_arkimex_type myark2 -ts_adapt_type none
268c4762a1bSJed Brown       requires: !single
269c4762a1bSJed Brown 
270c4762a1bSJed Brown TEST*/
271