xref: /petsc/src/ts/tutorials/ex16fwd.c (revision d5b43468fb8780a8feea140ccd6fa3e6a50411cc)
1 static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n\
2 Input parameters include:\n\
3       -mu : stiffness parameter\n\n";
4 
5 /* ------------------------------------------------------------------------
6 
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) = 0,
11    and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with an explicit Runge-Kutta method and its discrete tangent linear model.
12 
13    Notes:
14    This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t).
15 
16    (1) can be turned into a system of first order ODEs
17    [ y' ] = [          z          ]
18    [ z' ]   [ \mu (1 - y^2) z - y ]
19 
20    which then we can write as a vector equation
21 
22    [ u_1' ] = [             u_2           ]  (2)
23    [ u_2' ]   [ \mu (1 - u_1^2) u_2 - u_1 ]
24 
25    which is now in the form of u_t = F(u,t).
26 
27    The user provides the right-hand-side function
28 
29    [ f(u,t) ] = [ u_2                       ]
30                 [ \mu (1 - u_1^2) u_2 - u_1 ]
31 
32    the Jacobian function
33 
34    df   [       0           ;         1        ]
35    -- = [                                      ]
36    du   [ -2 \mu u_1*u_2 - 1;  \mu (1 - u_1^2) ]
37 
38    and the JacobainP (the Jacobian w.r.t. parameter) function
39 
40    df      [  0;   0;     0             ]
41    ---   = [                            ]
42    d\mu    [  0;   0;  (1 - u_1^2) u_2  ]
43 
44   ------------------------------------------------------------------------- */
45 
46 #include <petscts.h>
47 #include <petscmat.h>
48 typedef struct _n_User *User;
49 struct _n_User {
50   PetscReal mu;
51   PetscReal next_output;
52   PetscReal tprev;
53 };
54 
55 /*
56    User-defined routines
57 */
58 static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, void *ctx)
59 {
60   User               user = (User)ctx;
61   PetscScalar       *f;
62   const PetscScalar *x;
63 
64   PetscFunctionBeginUser;
65   PetscCall(VecGetArrayRead(X, &x));
66   PetscCall(VecGetArray(F, &f));
67   f[0] = x[1];
68   f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0];
69   PetscCall(VecRestoreArrayRead(X, &x));
70   PetscCall(VecRestoreArray(F, &f));
71   PetscFunctionReturn(0);
72 }
73 
74 static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, void *ctx)
75 {
76   User               user     = (User)ctx;
77   PetscReal          mu       = user->mu;
78   PetscInt           rowcol[] = {0, 1};
79   PetscScalar        J[2][2];
80   const PetscScalar *x;
81 
82   PetscFunctionBeginUser;
83   PetscCall(VecGetArrayRead(X, &x));
84   J[0][0] = 0;
85   J[1][0] = -2. * mu * x[1] * x[0] - 1.;
86   J[0][1] = 1.0;
87   J[1][1] = mu * (1.0 - x[0] * x[0]);
88   PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
89   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
90   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
91   if (A != B) {
92     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
93     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
94   }
95   PetscCall(VecRestoreArrayRead(X, &x));
96   PetscFunctionReturn(0);
97 }
98 
99 static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx)
100 {
101   PetscInt           row[] = {0, 1}, col[] = {2};
102   PetscScalar        J[2][1];
103   const PetscScalar *x;
104 
105   PetscFunctionBeginUser;
106   PetscCall(VecGetArrayRead(X, &x));
107   J[0][0] = 0;
108   J[1][0] = (1. - x[0] * x[0]) * x[1];
109   PetscCall(VecRestoreArrayRead(X, &x));
110   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
111 
112   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
113   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
114   PetscFunctionReturn(0);
115 }
116 
117 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
118 static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx)
119 {
120   const PetscScalar *x;
121   PetscReal          tfinal, dt, tprev;
122   User               user = (User)ctx;
123 
124   PetscFunctionBeginUser;
125   PetscCall(TSGetTimeStep(ts, &dt));
126   PetscCall(TSGetMaxTime(ts, &tfinal));
127   PetscCall(TSGetPrevTime(ts, &tprev));
128   PetscCall(VecGetArrayRead(X, &x));
129   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])));
130   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev));
131   PetscCall(VecRestoreArrayRead(X, &x));
132   PetscFunctionReturn(0);
133 }
134 
135 int main(int argc, char **argv)
136 {
137   TS             ts;   /* nonlinear solver */
138   Vec            x;    /* solution, residual vectors */
139   Mat            A;    /* Jacobian matrix */
140   Mat            Jacp; /* JacobianP matrix */
141   PetscInt       steps;
142   PetscReal      ftime   = 0.5;
143   PetscBool      monitor = PETSC_FALSE;
144   PetscScalar   *x_ptr;
145   PetscMPIInt    size;
146   struct _n_User user;
147   Mat            sp;
148 
149   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150      Initialize program
151      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152   PetscFunctionBeginUser;
153   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
154   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
155   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
156 
157   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158     Set runtime options
159     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160   user.mu          = 1;
161   user.next_output = 0.0;
162 
163   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
164   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
165 
166   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167     Create necessary matrix and vectors, solve same ODE on every process
168     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
170   PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
171   PetscCall(MatSetFromOptions(A));
172   PetscCall(MatSetUp(A));
173   PetscCall(MatCreateVecs(A, &x, NULL));
174 
175   PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
176   PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 3));
177   PetscCall(MatSetFromOptions(Jacp));
178   PetscCall(MatSetUp(Jacp));
179 
180   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 3, NULL, &sp));
181   PetscCall(MatZeroEntries(sp));
182   PetscCall(MatShift(sp, 1.0));
183 
184   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185      Create timestepping solver context
186      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
188   PetscCall(TSSetType(ts, TSRK));
189   PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
190   /*   Set RHS Jacobian for the adjoint integration */
191   PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &user));
192   PetscCall(TSSetMaxTime(ts, ftime));
193   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
194   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
195   PetscCall(TSForwardSetSensitivities(ts, 3, sp));
196   PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &user));
197 
198   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199      Set initial conditions
200    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201   PetscCall(VecGetArray(x, &x_ptr));
202 
203   x_ptr[0] = 2;
204   x_ptr[1] = 0.66666654321;
205   PetscCall(VecRestoreArray(x, &x_ptr));
206   PetscCall(TSSetTimeStep(ts, .001));
207 
208   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209      Set runtime options
210    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211   PetscCall(TSSetFromOptions(ts));
212 
213   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214      Solve nonlinear system
215      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216   PetscCall(TSSolve(ts, x));
217   PetscCall(TSGetSolveTime(ts, &ftime));
218   PetscCall(TSGetStepNumber(ts, &steps));
219   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime));
220   PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
221 
222   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n"));
223   PetscCall(MatView(sp, PETSC_VIEWER_STDOUT_WORLD));
224 
225   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226      Free work space.  All PETSc objects should be destroyed when they
227      are no longer needed.
228    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229   PetscCall(MatDestroy(&A));
230   PetscCall(MatDestroy(&Jacp));
231   PetscCall(VecDestroy(&x));
232   PetscCall(MatDestroy(&sp));
233   PetscCall(TSDestroy(&ts));
234   PetscCall(PetscFinalize());
235   return 0;
236 }
237 
238 /*TEST
239 
240     test:
241       args: -monitor 0 -ts_adapt_type none
242 
243 TEST*/
244