Actual source code: ex19.c


  2: static char help[] = "Solves the van der Pol DAE.\n\
  3: Input parameters include:\n";

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

  7:    This program solves the van der Pol DAE
  8:        y' = -z = f(y,z)        (1)
  9:        0  = y-(z^3/3 - z) = g(y,z)
 10:    on the domain 0 <= x <= 1, with the boundary conditions
 11:        y(0) = -2, y'(0) = -2.355301397608119909925287735864250951918
 12:    This is a nonlinear equation.

 14:    Notes:
 15:    This code demonstrates the TS solver interface with the Van der Pol DAE,
 16:    namely it is the case when there is no RHS (meaning the RHS == 0), and the
 17:    equations are converted to two variants of linear problems, u_t = f(u,t),
 18:    namely turning (1) into a vector equation in terms of u,

 20:    [     y' + z      ] = [ 0 ]
 21:    [ (z^3/3 - z) - y ]   [ 0 ]

 23:    which then we can write as a vector equation

 25:    [      u_1' + u_2       ] = [ 0 ]  (2)
 26:    [ (u_2^3/3 - u_2) - u_1 ]   [ 0 ]

 28:    which is now in the desired form of u_t = f(u,t). As this is a DAE, and
 29:    there is no u_2', there is no need for a split,

 31:    so

 33:    [ F(u',u,t) ] = [ u_1' ] + [         u_2           ]
 34:                    [  0   ]   [ (u_2^3/3 - u_2) - u_1 ]

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

 39:               dF   dF
 40:    J(F) = a * -- - --
 41:               du'  du

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

 45:    dF   [ 1 ; 0 ]
 46:    -- = [       ]
 47:    du'  [ 0 ; 0 ]

 49:    dF   [  0 ;      1     ]
 50:    -- = [                 ]
 51:    du   [ -1 ; 1 - u_2^2  ]

 53:    Hence,

 55:           [ a ;    -1     ]
 56:    J(F) = [               ]
 57:           [ 1 ; u_2^2 - 1 ]

 59:   ------------------------------------------------------------------------- */

 61: #include <petscts.h>

 63: typedef struct _n_User *User;
 64: struct _n_User {
 65:   PetscReal next_output;
 66: };

 68: /*
 69:    User-defined routines
 70: */

 72: static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx)
 73: {
 74:   PetscScalar       *f;
 75:   const PetscScalar *x, *xdot;

 78:   VecGetArrayRead(X, &x);
 79:   VecGetArrayRead(Xdot, &xdot);
 80:   VecGetArray(F, &f);
 81:   f[0] = xdot[0] + x[1];
 82:   f[1] = (x[1] * x[1] * x[1] / 3.0 - x[1]) - x[0];
 83:   VecRestoreArrayRead(X, &x);
 84:   VecRestoreArrayRead(Xdot, &xdot);
 85:   VecRestoreArray(F, &f);
 86:   return 0;
 87: }

 89: static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx)
 90: {
 91:   PetscInt           rowcol[] = {0, 1};
 92:   PetscScalar        J[2][2];
 93:   const PetscScalar *x;

 96:   VecGetArrayRead(X, &x);
 97:   J[0][0] = a;
 98:   J[0][1] = -1.;
 99:   J[1][0] = 1.;
100:   J[1][1] = -1. + x[1] * x[1];
101:   MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES);
102:   VecRestoreArrayRead(X, &x);

104:   MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY);
105:   MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY);
106:   if (A != B) {
107:     MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY);
108:     MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY);
109:   }
110:   return 0;
111: }

113: static PetscErrorCode RegisterMyARK2(void)
114: {
116:   {
117:     const PetscReal A[3][3] =
118:       {
119:         {0,                      0,    0},
120:         {0.41421356237309504880, 0,    0},
121:         {0.75,                   0.25, 0}
122:     },
123:                     At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL;
124:     TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL);
125:   }
126:   return 0;
127: }

129: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
130: static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx)
131: {
132:   const PetscScalar *x;
133:   PetscReal          tfinal, dt;
134:   User               user = (User)ctx;
135:   Vec                interpolatedX;

138:   TSGetTimeStep(ts, &dt);
139:   TSGetMaxTime(ts, &tfinal);

141:   while (user->next_output <= t && user->next_output <= tfinal) {
142:     VecDuplicate(X, &interpolatedX);
143:     TSInterpolate(ts, user->next_output, interpolatedX);
144:     VecGetArrayRead(interpolatedX, &x);
145:     PetscPrintf(PETSC_COMM_WORLD, "[%.1f] %3" 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]));
146:     VecRestoreArrayRead(interpolatedX, &x);
147:     VecDestroy(&interpolatedX);
148:     user->next_output += PetscRealConstant(0.1);
149:   }
150:   return 0;
151: }

153: int main(int argc, char **argv)
154: {
155:   TS             ts; /* nonlinear solver */
156:   Vec            x;  /* solution, residual vectors */
157:   Mat            A;  /* Jacobian matrix */
158:   PetscInt       steps;
159:   PetscReal      ftime   = 0.5;
160:   PetscBool      monitor = PETSC_FALSE;
161:   PetscScalar   *x_ptr;
162:   PetscMPIInt    size;
163:   struct _n_User user;

165:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
166:      Initialize program
167:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169:   PetscInitialize(&argc, &argv, NULL, help);
170:   MPI_Comm_size(PETSC_COMM_WORLD, &size);

173:   RegisterMyARK2();

175:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176:     Set runtime options
177:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */

179:   user.next_output = 0.0;
180:   PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL);

182:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
183:     Create necessary matrix and vectors, solve same ODE on every process
184:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
185:   MatCreate(PETSC_COMM_WORLD, &A);
186:   MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2);
187:   MatSetFromOptions(A);
188:   MatSetUp(A);
189:   MatCreateVecs(A, &x, NULL);

191:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192:      Create timestepping solver context
193:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
194:   TSCreate(PETSC_COMM_WORLD, &ts);
195:   TSSetType(ts, TSBEULER);
196:   TSSetIFunction(ts, NULL, IFunction, &user);
197:   TSSetIJacobian(ts, A, A, IJacobian, &user);
198:   TSSetMaxTime(ts, ftime);
199:   TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER);
200:   if (monitor) TSMonitorSet(ts, Monitor, &user, NULL);

202:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203:      Set initial conditions
204:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205:   VecGetArray(x, &x_ptr);
206:   x_ptr[0] = -2;
207:   x_ptr[1] = -2.355301397608119909925287735864250951918;
208:   VecRestoreArray(x, &x_ptr);
209:   TSSetTimeStep(ts, .001);

211:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212:      Set runtime options
213:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214:   TSSetFromOptions(ts);

216:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217:      Solve nonlinear system
218:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219:   TSSolve(ts, x);
220:   TSGetSolveTime(ts, &ftime);
221:   TSGetStepNumber(ts, &steps);
222:   PetscPrintf(PETSC_COMM_WORLD, "steps %3" PetscInt_FMT ", ftime %g\n", steps, (double)ftime);
223:   VecView(x, PETSC_VIEWER_STDOUT_WORLD);

225:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226:      Free work space.  All PETSc objects should be destroyed when they
227:      are no longer needed.
228:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229:   MatDestroy(&A);
230:   VecDestroy(&x);
231:   TSDestroy(&ts);

233:   PetscFinalize();
234:   return 0;
235: }

237: /*TEST

239:    test:
240:       requires: !single
241:       suffix: a
242:       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp
243:       output_file: output/ex19_pi42.out

245:    test:
246:       requires: !single
247:       suffix: b
248:       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_filter PI42
249:       output_file: output/ex19_pi42.out

251:    test:
252:       requires: !single
253:       suffix: c
254:       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_pid 0.4,0.2
255:       output_file: output/ex19_pi42.out

257:    test:
258:       requires: !single
259:       suffix: bdf_reject
260:       args: -ts_type bdf -ts_dt 0.5 -ts_max_steps 1 -ts_max_reject {{0 1 2}separate_output} -ts_error_if_step_fails false -ts_adapt_monitor

262: TEST*/