xref: /petsc/src/ts/tutorials/ex19.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Solves the van der Pol DAE.\n\
3c4762a1bSJed Brown Input parameters include:\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /* ------------------------------------------------------------------------
6c4762a1bSJed Brown 
7c4762a1bSJed Brown    This program solves the van der Pol DAE
8c4762a1bSJed Brown        y' = -z = f(y,z)        (1)
9c4762a1bSJed Brown        0  = y-(z^3/3 - z) = g(y,z)
10c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
11c4762a1bSJed Brown        y(0) = -2, y'(0) = -2.355301397608119909925287735864250951918
12c4762a1bSJed Brown    This is a nonlinear equation.
13c4762a1bSJed Brown 
14c4762a1bSJed Brown    Notes:
15c4762a1bSJed Brown    This code demonstrates the TS solver interface with the Van der Pol DAE,
16c4762a1bSJed Brown    namely it is the case when there is no RHS (meaning the RHS == 0), and the
17c4762a1bSJed Brown    equations are converted to two variants of linear problems, u_t = f(u,t),
18c4762a1bSJed Brown    namely turning (1) into a vector equation in terms of u,
19c4762a1bSJed Brown 
20c4762a1bSJed Brown    [     y' + z      ] = [ 0 ]
21c4762a1bSJed Brown    [ (z^3/3 - z) - y ]   [ 0 ]
22c4762a1bSJed Brown 
23c4762a1bSJed Brown    which then we can write as a vector equation
24c4762a1bSJed Brown 
25c4762a1bSJed Brown    [      u_1' + u_2       ] = [ 0 ]  (2)
26c4762a1bSJed Brown    [ (u_2^3/3 - u_2) - u_1 ]   [ 0 ]
27c4762a1bSJed Brown 
28c4762a1bSJed Brown    which is now in the desired form of u_t = f(u,t). As this is a DAE, and
29c4762a1bSJed Brown    there is no u_2', there is no need for a split,
30c4762a1bSJed Brown 
31c4762a1bSJed Brown    so
32c4762a1bSJed Brown 
335ab1ac2bSHong Zhang    [ F(u',u,t) ] = [ u_1' ] + [         u_2           ]
34c4762a1bSJed Brown                    [  0   ]   [ (u_2^3/3 - u_2) - u_1 ]
35c4762a1bSJed Brown 
365ab1ac2bSHong Zhang    Using the definition of the Jacobian of F (from the PETSc user manual),
375ab1ac2bSHong Zhang    in the equation F(u',u,t) = G(u,t),
38c4762a1bSJed Brown 
395ab1ac2bSHong Zhang               dF   dF
405ab1ac2bSHong Zhang    J(F) = a * -- - --
41c4762a1bSJed Brown               du'  du
42c4762a1bSJed Brown 
43c4762a1bSJed Brown    where d is the partial derivative. In this example,
44c4762a1bSJed Brown 
455ab1ac2bSHong Zhang    dF   [ 1 ; 0 ]
46c4762a1bSJed Brown    -- = [       ]
47c4762a1bSJed Brown    du'  [ 0 ; 0 ]
48c4762a1bSJed Brown 
495ab1ac2bSHong Zhang    dF   [  0 ;      1     ]
50c4762a1bSJed Brown    -- = [                 ]
51c4762a1bSJed Brown    du   [ -1 ; 1 - u_2^2  ]
52c4762a1bSJed Brown 
53c4762a1bSJed Brown    Hence,
54c4762a1bSJed Brown 
55c4762a1bSJed Brown           [ a ;    -1     ]
565ab1ac2bSHong Zhang    J(F) = [               ]
57c4762a1bSJed Brown           [ 1 ; u_2^2 - 1 ]
58c4762a1bSJed Brown 
59c4762a1bSJed Brown   ------------------------------------------------------------------------- */
60c4762a1bSJed Brown 
61c4762a1bSJed Brown #include <petscts.h>
62c4762a1bSJed Brown 
63c4762a1bSJed Brown typedef struct _n_User *User;
64c4762a1bSJed Brown struct _n_User {
65c4762a1bSJed Brown   PetscReal next_output;
66c4762a1bSJed Brown };
67c4762a1bSJed Brown 
68c4762a1bSJed Brown /*
690e3d61c9SBarry Smith    User-defined routines
70c4762a1bSJed Brown */
71c4762a1bSJed Brown 
729371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
73c4762a1bSJed Brown   PetscScalar       *f;
74c4762a1bSJed Brown   const PetscScalar *x, *xdot;
75c4762a1bSJed Brown 
76c4762a1bSJed Brown   PetscFunctionBeginUser;
779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
789566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
80c4762a1bSJed Brown   f[0] = xdot[0] + x[1];
81c4762a1bSJed Brown   f[1] = (x[1] * x[1] * x[1] / 3.0 - x[1]) - x[0];
829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
839566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
849566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
85c4762a1bSJed Brown   PetscFunctionReturn(0);
86c4762a1bSJed Brown }
87c4762a1bSJed Brown 
889371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
89c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
90c4762a1bSJed Brown   PetscScalar        J[2][2];
91c4762a1bSJed Brown   const PetscScalar *x;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown   PetscFunctionBeginUser;
949566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
959371c9d4SSatish Balay   J[0][0] = a;
969371c9d4SSatish Balay   J[0][1] = -1.;
979371c9d4SSatish Balay   J[1][0] = 1.;
989371c9d4SSatish Balay   J[1][1] = -1. + x[1] * x[1];
999566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
1009566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
101c4762a1bSJed Brown 
1029566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1039566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
104c4762a1bSJed Brown   if (A != B) {
1059566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
1069566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
107c4762a1bSJed Brown   }
108c4762a1bSJed Brown   PetscFunctionReturn(0);
109c4762a1bSJed Brown }
110c4762a1bSJed Brown 
1119371c9d4SSatish Balay static PetscErrorCode RegisterMyARK2(void) {
112c4762a1bSJed Brown   PetscFunctionBeginUser;
113c4762a1bSJed Brown   {
1149371c9d4SSatish Balay     const PetscReal A[3][3] =
1159371c9d4SSatish Balay       {
1169371c9d4SSatish Balay         {0,                      0,    0},
117c4762a1bSJed Brown         {0.41421356237309504880, 0,    0},
1189371c9d4SSatish Balay         {0.75,                   0.25, 0}
1199371c9d4SSatish Balay     },
1209371c9d4SSatish Balay                     At[3][3] = {{0, 0, 0}, {0.12132034355964257320, 0.29289321881345247560, 0}, {0.20710678118654752440, 0.50000000000000000000, 0.29289321881345247560}}, *bembedt = NULL, *bembed = NULL;
1219566063dSJacob Faibussowitsch     PetscCall(TSARKIMEXRegister("myark2", 2, 3, &At[0][0], NULL, NULL, &A[0][0], NULL, NULL, bembedt, bembed, 0, NULL, NULL));
122c4762a1bSJed Brown   }
123c4762a1bSJed Brown   PetscFunctionReturn(0);
124c4762a1bSJed Brown }
125c4762a1bSJed Brown 
126c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
1279371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) {
128c4762a1bSJed Brown   const PetscScalar *x;
129c4762a1bSJed Brown   PetscReal          tfinal, dt;
130c4762a1bSJed Brown   User               user = (User)ctx;
131c4762a1bSJed Brown   Vec                interpolatedX;
132c4762a1bSJed Brown 
133c4762a1bSJed Brown   PetscFunctionBeginUser;
1349566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1359566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
136c4762a1bSJed Brown 
137c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
1389566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &interpolatedX));
1399566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
1409566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX, &x));
14163a3b9bcSJacob Faibussowitsch     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])));
1429566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX, &x));
1439566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
144c4762a1bSJed Brown     user->next_output += PetscRealConstant(0.1);
145c4762a1bSJed Brown   }
146c4762a1bSJed Brown   PetscFunctionReturn(0);
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
1499371c9d4SSatish Balay int main(int argc, char **argv) {
150c4762a1bSJed Brown   TS             ts; /* nonlinear solver */
151c4762a1bSJed Brown   Vec            x;  /* solution, residual vectors */
152c4762a1bSJed Brown   Mat            A;  /* Jacobian matrix */
153c4762a1bSJed Brown   PetscInt       steps;
154c4762a1bSJed Brown   PetscReal      ftime   = 0.5;
155c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
156c4762a1bSJed Brown   PetscScalar   *x_ptr;
157c4762a1bSJed Brown   PetscMPIInt    size;
158c4762a1bSJed Brown   struct _n_User user;
159c4762a1bSJed Brown 
160c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161c4762a1bSJed Brown      Initialize program
162c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163327415f7SBarry Smith   PetscFunctionBeginUser;
1649566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1659566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1663c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
167c4762a1bSJed Brown 
1689566063dSJacob Faibussowitsch   PetscCall(RegisterMyARK2());
169c4762a1bSJed Brown 
170c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171c4762a1bSJed Brown     Set runtime options
172c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173c4762a1bSJed Brown 
174c4762a1bSJed Brown   user.next_output = 0.0;
1759566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
176c4762a1bSJed Brown 
177c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
179c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1809566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1819566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
1829566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1839566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
1849566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &x, NULL));
185c4762a1bSJed Brown 
186c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
187c4762a1bSJed Brown      Create timestepping solver context
188c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1899566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1909566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
1919566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
1929566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, IJacobian, &user));
1939566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, ftime));
1949566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
195*48a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
196c4762a1bSJed Brown 
197c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198c4762a1bSJed Brown      Set initial conditions
199c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2009566063dSJacob Faibussowitsch   PetscCall(VecGetArray(x, &x_ptr));
2019371c9d4SSatish Balay   x_ptr[0] = -2;
2029371c9d4SSatish Balay   x_ptr[1] = -2.355301397608119909925287735864250951918;
2039566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(x, &x_ptr));
2049566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, .001));
205c4762a1bSJed Brown 
206c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207c4762a1bSJed Brown      Set runtime options
208c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2099566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212c4762a1bSJed Brown      Solve nonlinear system
213c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2149566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, x));
2159566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2169566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
21763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "steps %3" PetscInt_FMT ", ftime %g\n", steps, (double)ftime));
2189566063dSJacob Faibussowitsch   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
219c4762a1bSJed Brown 
220c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
222c4762a1bSJed Brown      are no longer needed.
223c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&x));
2269566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
227c4762a1bSJed Brown 
2289566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
229b122ec5aSJacob Faibussowitsch   return 0;
230c4762a1bSJed Brown }
231c4762a1bSJed Brown 
232c4762a1bSJed Brown /*TEST
233c4762a1bSJed Brown 
234c4762a1bSJed Brown    test:
235c4762a1bSJed Brown       requires: !single
236c4762a1bSJed Brown       suffix: a
237c4762a1bSJed Brown       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp
238c4762a1bSJed Brown       output_file: output/ex19_pi42.out
239c4762a1bSJed Brown 
240c4762a1bSJed Brown    test:
241c4762a1bSJed Brown       requires: !single
242c4762a1bSJed Brown       suffix: b
243c4762a1bSJed Brown       args: -monitor -ts_type bdf -ts_rtol 1e-6 -ts_atol 1e-6 -ts_view -ts_adapt_type dsp -ts_adapt_dsp_filter PI42
244c4762a1bSJed Brown       output_file: output/ex19_pi42.out
245c4762a1bSJed Brown 
246c4762a1bSJed Brown    test:
247c4762a1bSJed Brown       requires: !single
248c4762a1bSJed Brown       suffix: c
249c4762a1bSJed Brown       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
250c4762a1bSJed Brown       output_file: output/ex19_pi42.out
251c4762a1bSJed Brown 
252e5b8ffdfSLisandro Dalcin    test:
253e5b8ffdfSLisandro Dalcin       requires: !single
254e5b8ffdfSLisandro Dalcin       suffix: bdf_reject
255e5b8ffdfSLisandro Dalcin       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
256e5b8ffdfSLisandro Dalcin 
257c4762a1bSJed Brown TEST*/
258