1c4762a1bSJed Brown static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n\
2c4762a1bSJed Brown Input parameters include:\n\
3c4762a1bSJed Brown -mu : stiffness parameter\n\n";
4c4762a1bSJed Brown
5c4762a1bSJed Brown /* ------------------------------------------------------------------------
6c4762a1bSJed Brown
7c4762a1bSJed Brown This program solves the van der Pol equation
8c4762a1bSJed Brown y'' - \mu (1-y^2)*y' + y = 0 (1)
9c4762a1bSJed Brown on the domain 0 <= x <= 1, with the boundary conditions
10c4762a1bSJed Brown y(0) = 2, y'(0) = 0,
11c4762a1bSJed Brown 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.
12c4762a1bSJed Brown
13c4762a1bSJed Brown Notes:
145ab1ac2bSHong Zhang This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t).
15c4762a1bSJed Brown
16c4762a1bSJed Brown (1) can be turned into a system of first order ODEs
17c4762a1bSJed Brown [ y' ] = [ z ]
18c4762a1bSJed Brown [ z' ] [ \mu (1 - y^2) z - y ]
19c4762a1bSJed Brown
20c4762a1bSJed Brown which then we can write as a vector equation
21c4762a1bSJed Brown
22c4762a1bSJed Brown [ u_1' ] = [ u_2 ] (2)
23c4762a1bSJed Brown [ u_2' ] [ \mu (1 - u_1^2) u_2 - u_1 ]
24c4762a1bSJed Brown
25c4762a1bSJed Brown which is now in the form of u_t = F(u,t).
26c4762a1bSJed Brown
27c4762a1bSJed Brown The user provides the right-hand-side function
28c4762a1bSJed Brown
295ab1ac2bSHong Zhang [ f(u,t) ] = [ u_2 ]
30c4762a1bSJed Brown [ \mu (1 - u_1^2) u_2 - u_1 ]
31c4762a1bSJed Brown
32c4762a1bSJed Brown the Jacobian function
33c4762a1bSJed Brown
345ab1ac2bSHong Zhang df [ 0 ; 1 ]
35c4762a1bSJed Brown -- = [ ]
36c4762a1bSJed Brown du [ -2 \mu u_1*u_2 - 1; \mu (1 - u_1^2) ]
37c4762a1bSJed Brown
38c4762a1bSJed Brown and the JacobainP (the Jacobian w.r.t. parameter) function
39c4762a1bSJed Brown
405ab1ac2bSHong Zhang df [ 0; 0; 0 ]
41c4762a1bSJed Brown --- = [ ]
42c4762a1bSJed Brown d\mu [ 0; 0; (1 - u_1^2) u_2 ]
43c4762a1bSJed Brown
44c4762a1bSJed Brown ------------------------------------------------------------------------- */
45c4762a1bSJed Brown
46c4762a1bSJed Brown #include <petscts.h>
47c4762a1bSJed Brown #include <petscmat.h>
48c4762a1bSJed Brown typedef struct _n_User *User;
49c4762a1bSJed Brown struct _n_User {
50c4762a1bSJed Brown PetscReal mu;
51c4762a1bSJed Brown PetscReal next_output;
52c4762a1bSJed Brown PetscReal tprev;
53c4762a1bSJed Brown };
54c4762a1bSJed Brown
55c4762a1bSJed Brown /*
560e3d61c9SBarry Smith User-defined routines
57c4762a1bSJed Brown */
RHSFunction(TS ts,PetscReal t,Vec X,Vec F,PetscCtx ctx)58*2a8381b2SBarry Smith static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec X, Vec F, PetscCtx ctx)
59d71ae5a4SJacob Faibussowitsch {
60c4762a1bSJed Brown User user = (User)ctx;
61c4762a1bSJed Brown PetscScalar *f;
62c4762a1bSJed Brown const PetscScalar *x;
63c4762a1bSJed Brown
64c4762a1bSJed Brown PetscFunctionBeginUser;
659566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
669566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f));
67c4762a1bSJed Brown f[0] = x[1];
68c4762a1bSJed Brown f[1] = user->mu * (1. - x[0] * x[0]) * x[1] - x[0];
699566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
709566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f));
713ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
72c4762a1bSJed Brown }
73c4762a1bSJed Brown
RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,PetscCtx ctx)74*2a8381b2SBarry Smith static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec X, Mat A, Mat B, PetscCtx ctx)
75d71ae5a4SJacob Faibussowitsch {
76c4762a1bSJed Brown User user = (User)ctx;
77c4762a1bSJed Brown PetscReal mu = user->mu;
78c4762a1bSJed Brown PetscInt rowcol[] = {0, 1};
79c4762a1bSJed Brown PetscScalar J[2][2];
80c4762a1bSJed Brown const PetscScalar *x;
81c4762a1bSJed Brown
82c4762a1bSJed Brown PetscFunctionBeginUser;
839566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
84c4762a1bSJed Brown J[0][0] = 0;
85c4762a1bSJed Brown J[1][0] = -2. * mu * x[1] * x[0] - 1.;
86c4762a1bSJed Brown J[0][1] = 1.0;
87c4762a1bSJed Brown J[1][1] = mu * (1.0 - x[0] * x[0]);
889566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
899566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
909566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
91c4762a1bSJed Brown if (A != B) {
929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
94c4762a1bSJed Brown }
959566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
963ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
97c4762a1bSJed Brown }
98c4762a1bSJed Brown
RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,PetscCtx ctx)99*2a8381b2SBarry Smith static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, PetscCtx ctx)
100d71ae5a4SJacob Faibussowitsch {
101c4762a1bSJed Brown PetscInt row[] = {0, 1}, col[] = {2};
102c4762a1bSJed Brown PetscScalar J[2][1];
103c4762a1bSJed Brown const PetscScalar *x;
104c4762a1bSJed Brown
105c4762a1bSJed Brown PetscFunctionBeginUser;
1069566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
107c4762a1bSJed Brown J[0][0] = 0;
108c4762a1bSJed Brown J[1][0] = (1. - x[0] * x[0]) * x[1];
1099566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
1109566063dSJacob Faibussowitsch PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
111c4762a1bSJed Brown
1129566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
1139566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
1143ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
115c4762a1bSJed Brown }
116c4762a1bSJed Brown
117c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
Monitor(TS ts,PetscInt step,PetscReal t,Vec X,PetscCtx ctx)118*2a8381b2SBarry Smith static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, PetscCtx ctx)
119d71ae5a4SJacob Faibussowitsch {
120c4762a1bSJed Brown const PetscScalar *x;
121c4762a1bSJed Brown PetscReal tfinal, dt, tprev;
122c4762a1bSJed Brown User user = (User)ctx;
123c4762a1bSJed Brown
124c4762a1bSJed Brown PetscFunctionBeginUser;
1259566063dSJacob Faibussowitsch PetscCall(TSGetTimeStep(ts, &dt));
1269566063dSJacob Faibussowitsch PetscCall(TSGetMaxTime(ts, &tfinal));
1279566063dSJacob Faibussowitsch PetscCall(TSGetPrevTime(ts, &tprev));
1289566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x));
12963a3b9bcSJacob Faibussowitsch 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])));
1309566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "t %.6f (tprev = %.6f) \n", (double)t, (double)tprev));
1319566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x));
1323ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
133c4762a1bSJed Brown }
134c4762a1bSJed Brown
main(int argc,char ** argv)135d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
136d71ae5a4SJacob Faibussowitsch {
137c4762a1bSJed Brown TS ts; /* nonlinear solver */
138c4762a1bSJed Brown Vec x; /* solution, residual vectors */
139c4762a1bSJed Brown Mat A; /* Jacobian matrix */
140c4762a1bSJed Brown Mat Jacp; /* JacobianP matrix */
141c4762a1bSJed Brown PetscInt steps;
142c4762a1bSJed Brown PetscReal ftime = 0.5;
143c4762a1bSJed Brown PetscBool monitor = PETSC_FALSE;
144c4762a1bSJed Brown PetscScalar *x_ptr;
145c4762a1bSJed Brown PetscMPIInt size;
146c4762a1bSJed Brown struct _n_User user;
147c4762a1bSJed Brown Mat sp;
148c4762a1bSJed Brown
149c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150c4762a1bSJed Brown Initialize program
151c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152327415f7SBarry Smith PetscFunctionBeginUser;
1539566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help));
1549566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1553c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
156c4762a1bSJed Brown
157c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
158c4762a1bSJed Brown Set runtime options
159c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
160c4762a1bSJed Brown user.mu = 1;
161c4762a1bSJed Brown user.next_output = 0.0;
162c4762a1bSJed Brown
1639566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
1649566063dSJacob Faibussowitsch PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
165c4762a1bSJed Brown
166c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167c4762a1bSJed Brown Create necessary matrix and vectors, solve same ODE on every process
168c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1699566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1709566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
1719566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A));
1729566063dSJacob Faibussowitsch PetscCall(MatSetUp(A));
1739566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &x, NULL));
174c4762a1bSJed Brown
1759566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &Jacp));
1769566063dSJacob Faibussowitsch PetscCall(MatSetSizes(Jacp, PETSC_DECIDE, PETSC_DECIDE, 2, 3));
1779566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(Jacp));
1789566063dSJacob Faibussowitsch PetscCall(MatSetUp(Jacp));
179c4762a1bSJed Brown
1809566063dSJacob Faibussowitsch PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, 2, 3, NULL, &sp));
1819566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(sp));
1829566063dSJacob Faibussowitsch PetscCall(MatShift(sp, 1.0));
183c4762a1bSJed Brown
184c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185c4762a1bSJed Brown Create timestepping solver context
186c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1879566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1889566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSRK));
1899566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, RHSFunction, &user));
190c4762a1bSJed Brown /* Set RHS Jacobian for the adjoint integration */
1919566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, RHSJacobian, &user));
1929566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, ftime));
1939566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
19448a46eb9SPierre Jolivet if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
1959566063dSJacob Faibussowitsch PetscCall(TSForwardSetSensitivities(ts, 3, sp));
1969566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobianP(ts, Jacp, RHSJacobianP, &user));
197c4762a1bSJed Brown
198c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199c4762a1bSJed Brown Set initial conditions
200c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2019566063dSJacob Faibussowitsch PetscCall(VecGetArray(x, &x_ptr));
202c4762a1bSJed Brown
2039371c9d4SSatish Balay x_ptr[0] = 2;
2049371c9d4SSatish Balay x_ptr[1] = 0.66666654321;
2059566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(x, &x_ptr));
2069566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .001));
207c4762a1bSJed Brown
208c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209c4762a1bSJed Brown Set runtime options
210c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2119566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts));
212c4762a1bSJed Brown
213c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214c4762a1bSJed Brown Solve nonlinear system
215c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2169566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, x));
2179566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime));
2189566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps));
21963a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, steps, (double)ftime));
2209566063dSJacob Faibussowitsch PetscCall(VecView(x, PETSC_VIEWER_STDOUT_WORLD));
221c4762a1bSJed Brown
2229566063dSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n"));
2239566063dSJacob Faibussowitsch PetscCall(MatView(sp, PETSC_VIEWER_STDOUT_WORLD));
224c4762a1bSJed Brown
225c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they
227c4762a1bSJed Brown are no longer needed.
228c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2299566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A));
2309566063dSJacob Faibussowitsch PetscCall(MatDestroy(&Jacp));
2319566063dSJacob Faibussowitsch PetscCall(VecDestroy(&x));
2329566063dSJacob Faibussowitsch PetscCall(MatDestroy(&sp));
2339566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts));
2349566063dSJacob Faibussowitsch PetscCall(PetscFinalize());
235b122ec5aSJacob Faibussowitsch return 0;
236c4762a1bSJed Brown }
237c4762a1bSJed Brown
238c4762a1bSJed Brown /*TEST
239c4762a1bSJed Brown
240c4762a1bSJed Brown test:
241c4762a1bSJed Brown args: -monitor 0 -ts_adapt_type none
242c4762a1bSJed Brown
243c4762a1bSJed Brown TEST*/
244