xref: /petsc/src/ts/tutorials/ex19.c (revision 9371c9d470a9602b6d10a8bf50c9b2280a79e45a)
1 
2 static char help[] = "Solves the van der Pol DAE.\n\
3 Input parameters include:\n";
4 
5 /* ------------------------------------------------------------------------
6 
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.
13 
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,
19 
20    [     y' + z      ] = [ 0 ]
21    [ (z^3/3 - z) - y ]   [ 0 ]
22 
23    which then we can write as a vector equation
24 
25    [      u_1' + u_2       ] = [ 0 ]  (2)
26    [ (u_2^3/3 - u_2) - u_1 ]   [ 0 ]
27 
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,
30 
31    so
32 
33    [ F(u',u,t) ] = [ u_1' ] + [         u_2           ]
34                    [  0   ]   [ (u_2^3/3 - u_2) - u_1 ]
35 
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),
38 
39               dF   dF
40    J(F) = a * -- - --
41               du'  du
42 
43    where d is the partial derivative. In this example,
44 
45    dF   [ 1 ; 0 ]
46    -- = [       ]
47    du'  [ 0 ; 0 ]
48 
49    dF   [  0 ;      1     ]
50    -- = [                 ]
51    du   [ -1 ; 1 - u_2^2  ]
52 
53    Hence,
54 
55           [ a ;    -1     ]
56    J(F) = [               ]
57           [ 1 ; u_2^2 - 1 ]
58 
59   ------------------------------------------------------------------------- */
60 
61 #include <petscts.h>
62 
63 typedef struct _n_User *User;
64 struct _n_User {
65   PetscReal next_output;
66 };
67 
68 /*
69    User-defined routines
70 */
71 
72 static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
73   PetscScalar       *f;
74   const PetscScalar *x, *xdot;
75 
76   PetscFunctionBeginUser;
77   PetscCall(VecGetArrayRead(X, &x));
78   PetscCall(VecGetArrayRead(Xdot, &xdot));
79   PetscCall(VecGetArray(F, &f));
80   f[0] = xdot[0] + x[1];
81   f[1] = (x[1] * x[1] * x[1] / 3.0 - x[1]) - x[0];
82   PetscCall(VecRestoreArrayRead(X, &x));
83   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
84   PetscCall(VecRestoreArray(F, &f));
85   PetscFunctionReturn(0);
86 }
87 
88 static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
89   PetscInt           rowcol[] = {0, 1};
90   PetscScalar        J[2][2];
91   const PetscScalar *x;
92 
93   PetscFunctionBeginUser;
94   PetscCall(VecGetArrayRead(X, &x));
95   J[0][0] = a;
96   J[0][1] = -1.;
97   J[1][0] = 1.;
98   J[1][1] = -1. + x[1] * x[1];
99   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
100   PetscCall(VecRestoreArrayRead(X, &x));
101 
102   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
103   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
104   if (A != B) {
105     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
106     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
107   }
108   PetscFunctionReturn(0);
109 }
110 
111 static PetscErrorCode RegisterMyARK2(void) {
112   PetscFunctionBeginUser;
113   {
114     const PetscReal A[3][3] =
115       {
116         {0,                      0,    0},
117         {0.41421356237309504880, 0,    0},
118         {0.75,                   0.25, 0}
119     },
120                     At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL;
121     PetscCall(TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL));
122   }
123   PetscFunctionReturn(0);
124 }
125 
126 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
127 static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) {
128   const PetscScalar *x;
129   PetscReal          tfinal, dt;
130   User               user = (User)ctx;
131   Vec                interpolatedX;
132 
133   PetscFunctionBeginUser;
134   PetscCall(TSGetTimeStep(ts, &dt));
135   PetscCall(TSGetMaxTime(ts, &tfinal));
136 
137   while (user->next_output <= t && user->next_output <= tfinal) {
138     PetscCall(VecDuplicate(X, &interpolatedX));
139     PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
140     PetscCall(VecGetArrayRead(interpolatedX, &x));
141     PetscCall(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])));
142     PetscCall(VecRestoreArrayRead(interpolatedX, &x));
143     PetscCall(VecDestroy(&interpolatedX));
144     user->next_output += PetscRealConstant(0.1);
145   }
146   PetscFunctionReturn(0);
147 }
148 
149 int main(int argc, char **argv) {
150   TS             ts; /* nonlinear solver */
151   Vec            x;  /* solution, residual vectors */
152   Mat            A;  /* Jacobian matrix */
153   PetscInt       steps;
154   PetscReal      ftime   = 0.5;
155   PetscBool      monitor = PETSC_FALSE;
156   PetscScalar   *x_ptr;
157   PetscMPIInt    size;
158   struct _n_User user;
159 
160   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161      Initialize program
162      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163   PetscFunctionBeginUser;
164   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
165   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
166   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
167 
168   PetscCall(RegisterMyARK2());
169 
170   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171     Set runtime options
172     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173 
174   user.next_output = 0.0;
175   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
176 
177   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178     Create necessary matrix and vectors, solve same ODE on every process
179     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
181   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
182   PetscCall(MatSetFromOptions(A));
183   PetscCall(MatSetUp(A));
184   PetscCall(MatCreateVecs(A, &x, NULL));
185 
186   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187      Create timestepping solver context
188      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
189   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
190   PetscCall(TSSetType(ts, TSBEULER));
191   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
192   PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user));
193   PetscCall(TSSetMaxTime(ts, ftime));
194   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
195   if (monitor) { PetscCall(TSMonitorSet(ts, Monitor, &user, NULL)); }
196 
197   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198      Set initial conditions
199    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200   PetscCall(VecGetArray(x, &x_ptr));
201   x_ptr[0] = -2;
202   x_ptr[1] = -2.355301397608119909925287735864250951918;
203   PetscCall(VecRestoreArray(x, &x_ptr));
204   PetscCall(TSSetTimeStep(ts, .001));
205 
206   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207      Set runtime options
208    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209   PetscCall(TSSetFromOptions(ts));
210 
211   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212      Solve nonlinear system
213      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214   PetscCall(TSSolve(ts, x));
215   PetscCall(TSGetSolveTime(ts, &ftime));
216   PetscCall(TSGetStepNumber(ts, &steps));
217   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %3" PetscInt_FMT ", ftime %g\n", steps, (double)ftime));
218   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
219 
220   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221      Free work space.  All PETSc objects should be destroyed when they
222      are no longer needed.
223    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224   PetscCall(MatDestroy(&A));
225   PetscCall(VecDestroy(&x));
226   PetscCall(TSDestroy(&ts));
227 
228   PetscCall(PetscFinalize());
229   return 0;
230 }
231 
232 /*TEST
233 
234    test:
235       requires: !single
236       suffix: a
237       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp
238       output_file: output/ex19_pi42.out
239 
240    test:
241       requires: !single
242       suffix: b
243       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_filter PI42
244       output_file: output/ex19_pi42.out
245 
246    test:
247       requires: !single
248       suffix: c
249       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
250       output_file: output/ex19_pi42.out
251 
252    test:
253       requires: !single
254       suffix: bdf_reject
255       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
256 
257 TEST*/
258