Actual source code: ex16.c

  1: static char help[] = "Solves the van der Pol equation and demonstrate IMEX.\n\
  2: Input parameters include:\n\
  3:       -mu : stiffness parameter\n\n";

  5: /* ------------------------------------------------------------------------

  7:    This program solves the van der Pol equation
  8:        y'' - \mu ((1-y^2)*y' - y) = 0        (1)
  9:    on the domain 0 <= x <= 1, with the boundary conditions
 10:        y(0) = 2, y'(0) = - 2/3 +10/(81*\mu) - 292/(2187*\mu^2),
 11:    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.

 13:    Notes:
 14:    This code demonstrates the TS solver interface to two variants of
 15:    linear problems, u_t = f(u,t), namely turning (1) into a system of
 16:    first order differential equations,

 18:    [ y' ] = [          z            ]
 19:    [ z' ]   [ \mu ((1 - y^2) z - y) ]

 21:    which then we can write as a vector equation

 23:    [ u_1' ] = [             u_2           ]  (2)
 24:    [ u_2' ]   [ \mu (1 - u_1^2) u_2 - u_1 ]

 26:    which is now in the desired form of u_t = f(u,t). One way that we
 27:    can split f(u,t) in (2) is to split by component,

 29:    [ u_1' ] = [ u_2 ] + [            0                ]
 30:    [ u_2' ]   [  0  ]   [ \mu ((1 - u_1^2) u_2 - u_1) ]

 32:    where

 34:    [ G(u,t) ] = [ u_2 ]
 35:                 [  0  ]

 37:    and

 39:    [ F(u',u,t) ] = [ u_1' ] - [            0                ]
 40:                    [ u_2' ]   [ \mu ((1 - u_1^2) u_2 - u_1) ]

 42:    Using the definition of the Jacobian of F (from the PETSc user manual),
 43:    in the equation F(u',u,t) = G(u,t),

 45:               dF   dF
 46:    J(F) = a * -- - --
 47:               du'  du

 49:    where d is the partial derivative. In this example,

 51:    dF   [ 1 ; 0 ]
 52:    -- = [       ]
 53:    du'  [ 0 ; 1 ]

 55:    dF   [       0             ;         0        ]
 56:    -- = [                                        ]
 57:    du   [ -\mu (2*u_1*u_2 + 1);  \mu (1 - u_1^2) ]

 59:    Hence,

 61:           [      a             ;          0          ]
 62:    J(F) = [                                          ]
 63:           [ \mu (2*u_1*u_2 + 1); a - \mu (1 - u_1^2) ]

 65:   ------------------------------------------------------------------------- */

 67: #include <petscts.h>

 69: typedef struct _n_User *User;
 70: struct _n_User {
 71:   PetscReal mu;
 72:   PetscBool imex;
 73:   PetscReal next_output;
 74: };

 76: /*
 77:    User-defined routines
 78: */
 79: static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx)
 80: {
 81:   User               user = (User)ctx;
 82:   PetscScalar       *f;
 83:   const PetscScalar *x;

 85:   PetscFunctionBeginUser;
 86:   PetscCall(VecGetArrayRead(X, &x));
 87:   PetscCall(VecGetArray(F, &f));
 88:   f[0] = (user->imex ? x[1] : 0);
 89:   f[1] = 0.0;
 90:   PetscCall(VecRestoreArrayRead(X, &x));
 91:   PetscCall(VecRestoreArray(F, &f));
 92:   PetscFunctionReturn(PETSC_SUCCESS);
 93: }

 95: static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
 96: {
 97:   User               user = (User)ctx;
 98:   const PetscScalar *x, *xdot;
 99:   PetscScalar       *f;

101:   PetscFunctionBeginUser;
102:   PetscCall(VecGetArrayRead(X, &x));
103:   PetscCall(VecGetArrayRead(Xdot, &xdot));
104:   PetscCall(VecGetArray(F, &f));
105:   f[0] = xdot[0] + (user->imex ? 0 : x[1]);
106:   f[1] = xdot[1] - user->mu * ((1. - x[0] * x[0]) * x[1] - x[0]);
107:   PetscCall(VecRestoreArrayRead(X, &x));
108:   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
109:   PetscCall(VecRestoreArray(F, &f));
110:   PetscFunctionReturn(PETSC_SUCCESS);
111: }

113: static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx)
114: {
115:   User               user     = (User)ctx;
116:   PetscReal          mu       = user->mu;
117:   PetscInt           rowcol[] = {0, 1};
118:   const PetscScalar *x;
119:   PetscScalar        J[2][2];

121:   PetscFunctionBeginUser;
122:   PetscCall(VecGetArrayRead(X, &x));
123:   J[0][0] = a;
124:   J[0][1] = (user->imex ? 0 : 1.);
125:   J[1][0] = mu * (2. * x[0] * x[1] + 1.);
126:   J[1][1] = a - mu * (1. - x[0] * x[0]);
127:   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
128:   PetscCall(VecRestoreArrayRead(X, &x));

130:   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
131:   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
132:   if (A != B) {
133:     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
134:     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
135:   }
136:   PetscFunctionReturn(PETSC_SUCCESS);
137: }

139: static PetscErrorCode RegisterMyARK2(void)
140: {
141:   PetscFunctionBeginUser;
142:   {
143:     const PetscReal A[3][3] =
144:       {
145:         {0,                      0,    0},
146:         {0.41421356237309504880, 0,    0},
147:         {0.75,                   0.25, 0}
148:     },
149:                     At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL;
150:     PetscCall(TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL));
151:   }
152:   PetscFunctionReturn(PETSC_SUCCESS);
153: }

155: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
156: static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx)
157: {
158:   const PetscScalar *x;
159:   PetscReal          tfinal, dt;
160:   User               user = (User)ctx;
161:   Vec                interpolatedX;

163:   PetscFunctionBeginUser;
164:   PetscCall(TSGetTimeStep(ts, &dt));
165:   PetscCall(TSGetMaxTime(ts, &tfinal));

167:   while (user->next_output <= t && user->next_output <= tfinal) {
168:     PetscCall(VecDuplicate(X, &interpolatedX));
169:     PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
170:     PetscCall(VecGetArrayRead(interpolatedX, &x));
171:     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])));
172:     PetscCall(VecRestoreArrayRead(interpolatedX, &x));
173:     PetscCall(VecDestroy(&interpolatedX));

175:     user->next_output += 0.1;
176:   }
177:   PetscFunctionReturn(PETSC_SUCCESS);
178: }

180: int main(int argc, char **argv)
181: {
182:   TS             ts; /* nonlinear solver */
183:   Vec            x;  /* solution, residual vectors */
184:   Mat            A;  /* Jacobian matrix */
185:   PetscInt       steps;
186:   PetscReal      ftime   = 0.5;
187:   PetscBool      monitor = PETSC_FALSE;
188:   PetscScalar   *x_ptr;
189:   PetscMPIInt    size;
190:   struct _n_User user;

192:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193:      Initialize program
194:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195:   PetscFunctionBeginUser;
196:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
197:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
198:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");

200:   PetscCall(RegisterMyARK2());

202:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203:     Set runtime options
204:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205:   user.mu          = 1000.0;
206:   user.imex        = PETSC_TRUE;
207:   user.next_output = 0.0;

209:   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
210:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-imex", &user.imex, NULL));
211:   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));

213:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214:     Create necessary matrix and vectors, solve same ODE on every process
215:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
217:   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
218:   PetscCall(MatSetFromOptions(A));
219:   PetscCall(MatSetUp(A));
220:   PetscCall(MatCreateVecs(A, &x, NULL));

222:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
223:      Create timestepping solver context
224:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
225:   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
226:   PetscCall(TSSetType(ts, TSBEULER));
227:   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
228:   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
229:   PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user));
230:   PetscCall(TSSetMaxTime(ts, ftime));
231:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
232:   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));

234:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
235:      Set initial conditions
236:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
237:   PetscCall(VecGetArray(x, &x_ptr));
238:   x_ptr[0] = 2.0;
239:   x_ptr[1] = -2.0 / 3.0 + 10.0 / (81.0 * user.mu) - 292.0 / (2187.0 * user.mu * user.mu);
240:   PetscCall(VecRestoreArray(x, &x_ptr));
241:   PetscCall(TSSetTimeStep(ts, 0.01));

243:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244:      Set runtime options
245:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246:   PetscCall(TSSetFromOptions(ts));

248:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
249:      Solve nonlinear system
250:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
251:   PetscCall(TSSolve(ts, x));
252:   PetscCall(TSGetSolveTime(ts, &ftime));
253:   PetscCall(TSGetStepNumber(ts, &steps));
254:   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime));
255:   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));

257:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
258:      Free work space.  All PETSc objects should be destroyed when they
259:      are no longer needed.
260:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
261:   PetscCall(MatDestroy(&A));
262:   PetscCall(VecDestroy(&x));
263:   PetscCall(TSDestroy(&ts));

265:   PetscCall(PetscFinalize());
266:   return 0;
267: }

269: /*TEST

271:     test:
272:       args: -ts_type arkimex -ts_arkimex_type myark2 -ts_adapt_type none
273:       requires: !single

275: TEST*/