1c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation and demonstrate IMEX.\n\
2c4762a1bSJed Brown Input parameters include:\n\
3c4762a1bSJed Brown -mu : stiffness parameter\n\n";
4c4762a1bSJed Brown
5c4762a1bSJed Brown /* ------------------------------------------------------------------------
6c4762a1bSJed Brown
7c4762a1bSJed Brown This program solves the van der Pol equation
8c4762a1bSJed Brown y'' - \mu ((1-y^2)*y' - y) = 0 (1)
9c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions
10c4762a1bSJed Brown y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2),
11c4762a1bSJed 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.
12c4762a1bSJed Brown
13c4762a1bSJed Brown Notes:
14c4762a1bSJed Brown This code demonstrates the TS solver interface to two variants of
15c4762a1bSJed Brown linear problems, u_t = f(u,t), namely turning (1) into a system of
16c4762a1bSJed Brown first order differential equations,
17c4762a1bSJed Brown
18c4762a1bSJed Brown [ y' ] = [ z ]
19c4762a1bSJed Brown [ z' ] [ \mu ((1 - y^2) z - y) ]
20c4762a1bSJed Brown
21c4762a1bSJed Brown which then we can write as a vector equation
22c4762a1bSJed Brown
23c4762a1bSJed Brown [ u_1' ] = [ u_2 ] (2)
24c4762a1bSJed Brown [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ]
25c4762a1bSJed Brown
26c4762a1bSJed Brown which is now in the desired form of u_t = f(u,t). One way that we
27c4762a1bSJed Brown can split f(u,t) in (2) is to split by component,
28c4762a1bSJed Brown
29c4762a1bSJed Brown [ u_1' ] = [ u_2 ] + [ 0 ]
30c4762a1bSJed Brown [ u_2' ] [ 0 ] [ \mu ((1 - u_1^2) u_2 - u_1) ]
31c4762a1bSJed Brown
32c4762a1bSJed Brown where
33c4762a1bSJed Brown
345ab1ac2bSHong Zhang [ G(u,t) ] = [ u_2 ]
35c4762a1bSJed Brown [ 0 ]
36c4762a1bSJed Brown
37c4762a1bSJed Brown and
38c4762a1bSJed Brown
395ab1ac2bSHong Zhang [ F(u',u,t) ] = [ u_1' ] - [ 0 ]
40c4762a1bSJed Brown [ u_2' ] [ \mu ((1 - u_1^2) u_2 - u_1) ]
41c4762a1bSJed Brown
425ab1ac2bSHong Zhang Using the definition of the Jacobian of F (from the PETSc user manual),
435ab1ac2bSHong Zhang in the equation F(u',u,t) = G(u,t),
44c4762a1bSJed Brown
455ab1ac2bSHong Zhang dF dF
465ab1ac2bSHong Zhang J(F) = a * -- - --
47c4762a1bSJed Brown du' du
48c4762a1bSJed Brown
49c4762a1bSJed Brown where d is the partial derivative. In this example,
50c4762a1bSJed Brown
515ab1ac2bSHong Zhang dF [ 1 ; 0 ]
52c4762a1bSJed Brown -- = [ ]
53c4762a1bSJed Brown du' [ 0 ; 1 ]
54c4762a1bSJed Brown
555ab1ac2bSHong Zhang dF [ 0 ; 0 ]
56c4762a1bSJed Brown -- = [ ]
57c4762a1bSJed Brown du [ -\mu (2*u_1*u_2 + 1); \mu (1 - u_1^2) ]
58c4762a1bSJed Brown
59c4762a1bSJed Brown Hence,
60c4762a1bSJed Brown
61c4762a1bSJed Brown [ a ; 0 ]
625ab1ac2bSHong Zhang J(F) = [ ]
63c4762a1bSJed Brown [ \mu (2*u_1*u_2 + 1); a - \mu (1 - u_1^2) ]
64c4762a1bSJed Brown
65c4762a1bSJed Brown ------------------------------------------------------------------------- */
66c4762a1bSJed Brown
67c4762a1bSJed Brown #include <petscts.h>
68c4762a1bSJed Brown
69c4762a1bSJed Brown typedef struct _n_User *User;
70c4762a1bSJed Brown struct _n_User {
71c4762a1bSJed Brown PetscReal mu;
72c4762a1bSJed Brown PetscBool imex;
73c4762a1bSJed Brown PetscReal next_output;
74c4762a1bSJed Brown };
75c4762a1bSJed Brown
76c4762a1bSJed Brown /*
770e3d61c9SBarry Smith User-defined routines
78c4762a1bSJed Brown */
RHSFunction(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)79*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
80d71ae5a4SJacob Faibussowitsch {
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));
923ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
93c4762a1bSJed Brown }
94c4762a1bSJed Brown
IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,PetscCtx ctx)95*2a8381b2SBarry Smith static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, PetscCtx ctx)
96d71ae5a4SJacob Faibussowitsch {
97c4762a1bSJed Brown User user = (User)ctx;
98c4762a1bSJed Brown const PetscScalar *x, *xdot;
99c4762a1bSJed Brown PetscScalar *f;
100c4762a1bSJed Brown
101c4762a1bSJed Brown PetscFunctionBeginUser;
1029566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
1039566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Xdot, &xdot));
1049566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
105c4762a1bSJed Brown f[0] = xdot[0] + (user->imex ? 0 : x[1]);
106c4762a1bSJed Brown f[1] = xdot[1] - user->mu * ((1. - x[0] * x[0]) * x[1] - x[0]);
1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
1089566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Xdot, &xdot));
1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
1103ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
111c4762a1bSJed Brown }
112c4762a1bSJed Brown
IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,PetscCtx ctx)113*2a8381b2SBarry Smith static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, PetscCtx ctx)
114d71ae5a4SJacob Faibussowitsch {
115c4762a1bSJed Brown User user = (User)ctx;
116c4762a1bSJed Brown PetscReal mu = user->mu;
117c4762a1bSJed Brown PetscInt rowcol[] = {0, 1};
118c4762a1bSJed Brown const PetscScalar *x;
119c4762a1bSJed Brown PetscScalar J[2][2];
120c4762a1bSJed Brown
121c4762a1bSJed Brown PetscFunctionBeginUser;
1229566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
1239371c9d4SSatish Balay J[0][0] = a;
1249371c9d4SSatish Balay J[0][1] = (user->imex ? 0 : 1.);
1259371c9d4SSatish Balay J[1][0] = mu * (2. * x[0] * x[1] + 1.);
1269371c9d4SSatish Balay J[1][1] = a - mu * (1. - x[0] * x[0]);
1279566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
1289566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
129c4762a1bSJed Brown
1309566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1319566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
132c4762a1bSJed Brown if (A != B) {
1339566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1349566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
135c4762a1bSJed Brown }
1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
137c4762a1bSJed Brown }
138c4762a1bSJed Brown
RegisterMyARK2(void)139d71ae5a4SJacob Faibussowitsch static PetscErrorCode RegisterMyARK2(void)
140d71ae5a4SJacob Faibussowitsch {
141c4762a1bSJed Brown PetscFunctionBeginUser;
142c4762a1bSJed Brown {
1439371c9d4SSatish Balay const PetscReal A[3][3] =
1449371c9d4SSatish Balay {
1459371c9d4SSatish Balay {0, 0, 0},
146c4762a1bSJed Brown {0.41421356237309504880, 0, 0},
1479371c9d4SSatish Balay {0.75, 0.25, 0}
1489371c9d4SSatish Balay },
1499371c9d4SSatish Balay At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL;
1509566063dSJacob Faibussowitsch PetscCall(TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL));
151c4762a1bSJed Brown }
1523ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
153c4762a1bSJed Brown }
154c4762a1bSJed Brown
155c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
Monitor(TS ts,PetscInt step,PetscReal t,Vec X,PetscCtx ctx)156*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, PetscCtx ctx)
157d71ae5a4SJacob Faibussowitsch {
158c4762a1bSJed Brown const PetscScalar *x;
159c4762a1bSJed Brown PetscReal tfinal, dt;
160c4762a1bSJed Brown User user = (User)ctx;
161c4762a1bSJed Brown Vec interpolatedX;
162c4762a1bSJed Brown
163c4762a1bSJed Brown PetscFunctionBeginUser;
1649566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt));
1659566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal));
166c4762a1bSJed Brown
167c4762a1bSJed Brown while (user->next_output <= t && user->next_output <= tfinal) {
1689566063dSJacob Faibussowitsch PetscCall(VecDuplicate(X, &interpolatedX));
1699566063dSJacob Faibussowitsch PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
1709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(interpolatedX, &x));
17163a3b9bcSJacob 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])));
1729566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(interpolatedX, &x));
1739566063dSJacob Faibussowitsch PetscCall(VecDestroy(&interpolatedX));
174c4762a1bSJed Brown
175c4762a1bSJed Brown user->next_output += 0.1;
176c4762a1bSJed Brown }
1773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
178c4762a1bSJed Brown }
179c4762a1bSJed Brown
main(int argc,char ** argv)180d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
181d71ae5a4SJacob Faibussowitsch {
182c4762a1bSJed Brown TS ts; /* nonlinear solver */
183c4762a1bSJed Brown Vec x; /* solution, residual vectors */
184c4762a1bSJed Brown Mat A; /* Jacobian matrix */
185c4762a1bSJed Brown PetscInt steps;
186c4762a1bSJed Brown PetscReal ftime = 0.5;
187c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE;
188c4762a1bSJed Brown PetscScalar *x_ptr;
189c4762a1bSJed Brown PetscMPIInt size;
190c4762a1bSJed Brown struct _n_User user;
191c4762a1bSJed Brown
192c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193c4762a1bSJed Brown Initialize program
194c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195327415f7SBarry Smith PetscFunctionBeginUser;
1969566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1979566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1983c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
199c4762a1bSJed Brown
2009566063dSJacob Faibussowitsch PetscCall(RegisterMyARK2());
201c4762a1bSJed Brown
202c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203c4762a1bSJed Brown Set runtime options
204c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205c4762a1bSJed Brown user.mu = 1000.0;
206c4762a1bSJed Brown user.imex = PETSC_TRUE;
207c4762a1bSJed Brown user.next_output = 0.0;
208c4762a1bSJed Brown
2099566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
2109566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-imex", &user.imex, NULL));
2119566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
212c4762a1bSJed Brown
213c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process
215c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2169566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
2179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
2189566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
2199566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
2209566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL));
221c4762a1bSJed Brown
222c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223c4762a1bSJed Brown Create timestepping solver context
224c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2259566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
2269566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSBEULER));
2279566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
2289566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
2299566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user));
2309566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime));
2319566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
23248a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
233c4762a1bSJed Brown
234c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235c4762a1bSJed Brown Set initial conditions
236c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2379566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr));
238c4762a1bSJed Brown x_ptr[0] = 2.0;
239c4762a1bSJed Brown x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
2409566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr));
2419566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.01));
242c4762a1bSJed Brown
243c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244c4762a1bSJed Brown Set runtime options
245c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2469566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
247c4762a1bSJed Brown
248c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249c4762a1bSJed Brown Solve nonlinear system
250c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2519566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
2529566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime));
2539566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
25463a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime));
2559566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
256c4762a1bSJed Brown
257c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
259c4762a1bSJed Brown are no longer needed.
260c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2619566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2629566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
2639566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
264c4762a1bSJed Brown
2659566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
266b122ec5aSJacob Faibussowitsch return 0;
267c4762a1bSJed Brown }
268c4762a1bSJed Brown
269c4762a1bSJed Brown /*TEST
270c4762a1bSJed Brown
271c4762a1bSJed Brown test:
272c4762a1bSJed Brown args: -ts_type arkimex -ts_arkimex_type myark2 -ts_adapt_type none
273c4762a1bSJed Brown requires: !single
274c4762a1bSJed Brown
275c4762a1bSJed Brown TEST*/
276