Actual source code: ex16.c
2: static char help[] = "Solves the van der Pol equation and demonstrate IMEX.\n\
3: Input parameters include:\n\
4: -mu : stiffness parameter\n\n";
6: /* ------------------------------------------------------------------------
8: This program solves the van der Pol equation
9: y'' - \mu ((1-y^2)*y' - y) = 0 (1)
10: on the domain 0 <= x <= 1, with the boundary conditions
11: y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2),
12: 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.
14: Notes:
15: This code demonstrates the TS solver interface to two variants of
16: linear problems, u_t = f(u,t), namely turning (1) into a system of
17: first order differential equations,
19: [ y' ] = [ z ]
20: [ z' ] [ \mu ((1 - y^2) z - y) ]
22: which then we can write as a vector equation
24: [ u_1' ] = [ u_2 ] (2)
25: [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ]
27: which is now in the desired form of u_t = f(u,t). One way that we
28: can split f(u,t) in (2) is to split by component,
30: [ u_1' ] = [ u_2 ] + [ 0 ]
31: [ u_2' ] [ 0 ] [ \mu ((1 - u_1^2) u_2 - u_1) ]
33: where
35: [ G(u,t) ] = [ u_2 ]
36: [ 0 ]
38: and
40: [ F(u',u,t) ] = [ u_1' ] - [ 0 ]
41: [ u_2' ] [ \mu ((1 - u_1^2) u_2 - u_1) ]
43: Using the definition of the Jacobian of F (from the PETSc user manual),
44: in the equation F(u',u,t) = G(u,t),
46: dF dF
47: J(F) = a * -- - --
48: du' du
50: where d is the partial derivative. In this example,
52: dF [ 1 ; 0 ]
53: -- = [ ]
54: du' [ 0 ; 1 ]
56: dF [ 0 ; 0 ]
57: -- = [ ]
58: du [ -\mu (2*u_1*u_2 + 1); \mu (1 - u_1^2) ]
60: Hence,
62: [ a ; 0 ]
63: J(F) = [ ]
64: [ \mu (2*u_1*u_2 + 1); a - \mu (1 - u_1^2) ]
66: ------------------------------------------------------------------------- */
68: #include <petscts.h>
70: typedef struct _n_User *User;
71: struct _n_User {
72: PetscReal mu;
73: PetscBool imex;
74: PetscReal next_output;
75: };
77: /*
78: User-defined routines
79: */
80: static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx)
81: {
82: User user = (User)ctx;
83: PetscScalar *f;
84: const PetscScalar *x;
87: VecGetArrayRead(X, &x);
88: VecGetArray(F, &f);
89: f[0] = (user->imex ? x[1] : 0);
90: f[1] = 0.0;
91: VecRestoreArrayRead(X, &x);
92: VecRestoreArray(F, &f);
93: return 0;
94: }
96: static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
97: {
98: User user = (User)ctx;
99: const PetscScalar *x, *xdot;
100: PetscScalar *f;
103: VecGetArrayRead(X, &x);
104: VecGetArrayRead(Xdot, &xdot);
105: VecGetArray(F, &f);
106: f[0] = xdot[0] + (user->imex ? 0 : x[1]);
107: f[1] = xdot[1] - user->mu * ((1. - x[0] * x[0]) * x[1] - x[0]);
108: VecRestoreArrayRead(X, &x);
109: VecRestoreArrayRead(Xdot, &xdot);
110: VecRestoreArray(F, &f);
111: return 0;
112: }
114: static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx)
115: {
116: User user = (User)ctx;
117: PetscReal mu = user->mu;
118: PetscInt rowcol[] = {0, 1};
119: const PetscScalar *x;
120: PetscScalar J[2][2];
123: VecGetArrayRead(X, &x);
124: J[0][0] = a;
125: J[0][1] = (user->imex ? 0 : 1.);
126: J[1][0] = mu * (2. * x[0] * x[1] + 1.);
127: J[1][1] = a - mu * (1. - x[0] * x[0]);
128: MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES);
129: VecRestoreArrayRead(X, &x);
131: MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
132: MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
133: if (A != B) {
134: MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
135: MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
136: }
137: return 0;
138: }
140: static PetscErrorCode RegisterMyARK2(void)
141: {
143: {
144: const PetscReal A[3][3] =
145: {
146: {0, 0, 0},
147: {0.41421356237309504880, 0, 0},
148: {0.75, 0.25, 0}
149: },
150: At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL;
151: TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL);
152: }
153: return 0;
154: }
156: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
157: static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx)
158: {
159: const PetscScalar *x;
160: PetscReal tfinal, dt;
161: User user = (User)ctx;
162: Vec interpolatedX;
165: TSGetTimeStep(ts, &dt);
166: TSGetMaxTime(ts, &tfinal);
168: while (user->next_output <= t && user->next_output <= tfinal) {
169: VecDuplicate(X, &interpolatedX);
170: TSInterpolate(ts, user->next_output, interpolatedX);
171: VecGetArrayRead(interpolatedX, &x);
172: 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]));
173: VecRestoreArrayRead(interpolatedX, &x);
174: VecDestroy(&interpolatedX);
176: user->next_output += 0.1;
177: }
178: return 0;
179: }
181: int main(int argc, char **argv)
182: {
183: TS ts; /* nonlinear solver */
184: Vec x; /* solution, residual vectors */
185: Mat A; /* Jacobian matrix */
186: PetscInt steps;
187: PetscReal ftime = 0.5;
188: PetscBool monitor = PETSC_FALSE;
189: PetscScalar *x_ptr;
190: PetscMPIInt size;
191: struct _n_User user;
193: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194: Initialize program
195: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197: PetscInitialize(&argc, &argv, NULL, help);
198: MPI_Comm_size(PETSC_COMM_WORLD, &size);
201: RegisterMyARK2();
203: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204: Set runtime options
205: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206: user.mu = 1000.0;
207: user.imex = PETSC_TRUE;
208: user.next_output = 0.0;
210: PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL);
211: PetscOptionsGetBool(NULL, NULL, "-imex", &user.imex, NULL);
212: PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL);
214: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215: Create necessary matrix and vectors, solve same ODE on every process
216: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217: MatCreate(PETSC_COMM_WORLD, &A);
218: MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2);
219: MatSetFromOptions(A);
220: MatSetUp(A);
221: MatCreateVecs(A, &x, NULL);
223: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224: Create timestepping solver context
225: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226: TSCreate(PETSC_COMM_WORLD, &ts);
227: TSSetType(ts, TSBEULER);
228: TSSetRHSFunction(ts, NULL, RHSFunction, &user);
229: TSSetIFunction(ts, NULL, IFunction, &user);
230: TSSetIJacobian(ts, A, A, IJacobian, &user);
231: TSSetMaxTime(ts, ftime);
232: TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);
233: if (monitor) TSMonitorSet(ts, Monitor, &user, NULL);
235: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
236: Set initial conditions
237: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
238: VecGetArray(x, &x_ptr);
239: x_ptr[0] = 2.0;
240: x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
241: VecRestoreArray(x, &x_ptr);
242: TSSetTimeStep(ts, 0.01);
244: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
245: Set runtime options
246: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
247: TSSetFromOptions(ts);
249: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
250: Solve nonlinear system
251: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
252: TSSolve(ts, x);
253: TSGetSolveTime(ts, &ftime);
254: TSGetStepNumber(ts, &steps);
255: PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime);
256: VecView(x, PETSC_VIEWER_STDOUT_WORLD);
258: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
259: Free work space. All PETSc objects should be destroyed when they
260: are no longer needed.
261: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
262: MatDestroy(&A);
263: VecDestroy(&x);
264: TSDestroy(&ts);
266: PetscFinalize();
267: return 0;
268: }
270: /*TEST
272: test:
273: args: -ts_type arkimex -ts_arkimex_type myark2 -ts_adapt_type none
274: requires: !single
276: TEST*/