xref: /petsc/src/ts/tutorials/ex20fwd.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
1c4762a1bSJed Brown #define c11 1.0
2c4762a1bSJed Brown #define c12 0
3c4762a1bSJed Brown #define c21 2.0
4c4762a1bSJed Brown #define c22 1.0
5c4762a1bSJed Brown static char help[] = "Solves the van der Pol equation.\n\
6c4762a1bSJed Brown Input parameters include:\n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown /* ------------------------------------------------------------------------
9c4762a1bSJed Brown 
10c4762a1bSJed Brown    This code demonstrates how to compute trajectory sensitivties w.r.t. the stiffness parameter mu.
112d4ee042Sprj-    1) Use two vectors s and sp for sensitivities w.r.t. initial values and paraeters respectively. This case requires the original Jacobian matrix and a JacobianP matrix for the only parameter mu.
122d4ee042Sprj-    2) Consider the initial values to be parameters as well. Then there are three parameters in total. The JacobianP matrix will be combined matrix of the Jacobian matrix and JacobianP matrix in the previous case. This choice can be selected by using command line option '-combined'
13c4762a1bSJed Brown 
14c4762a1bSJed Brown   ------------------------------------------------------------------------- */
15c4762a1bSJed Brown #include <petscts.h>
16c4762a1bSJed Brown #include <petsctao.h>
17c4762a1bSJed Brown 
18c4762a1bSJed Brown typedef struct _n_User *User;
19c4762a1bSJed Brown struct _n_User {
20c4762a1bSJed Brown   PetscReal mu;
21c4762a1bSJed Brown   PetscReal next_output;
22c4762a1bSJed Brown   PetscBool combined;
23c4762a1bSJed Brown   /* Sensitivity analysis support */
24c4762a1bSJed Brown   PetscInt  steps;
25c4762a1bSJed Brown   PetscReal ftime;
26c4762a1bSJed Brown   Mat       Jac;  /* Jacobian matrix */
27c4762a1bSJed Brown   Mat       Jacp; /* JacobianP matrix */
28c4762a1bSJed Brown   Vec       x;
29c4762a1bSJed Brown   Mat       sp; /* forward sensitivity variables */
30c4762a1bSJed Brown };
31c4762a1bSJed Brown 
32c4762a1bSJed Brown /*
330e3d61c9SBarry Smith    User-defined routines
34c4762a1bSJed Brown */
359371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) {
36c4762a1bSJed Brown   User               user = (User)ctx;
37c4762a1bSJed Brown   const PetscScalar *x, *xdot;
38c4762a1bSJed Brown   PetscScalar       *f;
39c4762a1bSJed Brown 
40c4762a1bSJed Brown   PetscFunctionBeginUser;
419566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
429566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot, &xdot));
439566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F, &f));
44c4762a1bSJed Brown   f[0] = xdot[0] - x[1];
45c4762a1bSJed Brown   f[1] = c21 * (xdot[0] - x[1]) + xdot[1] - user->mu * ((1.0 - x[0] * x[0]) * x[1] - x[0]);
469566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot, &xdot));
489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F, &f));
49c4762a1bSJed Brown   PetscFunctionReturn(0);
50c4762a1bSJed Brown }
51c4762a1bSJed Brown 
529371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, void *ctx) {
53c4762a1bSJed Brown   User               user     = (User)ctx;
54c4762a1bSJed Brown   PetscInt           rowcol[] = {0, 1};
55c4762a1bSJed Brown   PetscScalar        J[2][2];
56c4762a1bSJed Brown   const PetscScalar *x;
57c4762a1bSJed Brown 
58c4762a1bSJed Brown   PetscFunctionBeginUser;
599566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
609371c9d4SSatish Balay   J[0][0] = a;
619371c9d4SSatish Balay   J[0][1] = -1.0;
629371c9d4SSatish Balay   J[1][0] = c21 * a + user->mu * (1.0 + 2.0 * x[0] * x[1]);
639371c9d4SSatish Balay   J[1][1] = -c21 + a - user->mu * (1.0 - x[0] * x[0]);
649566063dSJacob Faibussowitsch   PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES));
659566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
66c4762a1bSJed Brown 
679566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
689566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
69c4762a1bSJed Brown   if (A != B) {
709566063dSJacob Faibussowitsch     PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY));
719566063dSJacob Faibussowitsch     PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY));
72c4762a1bSJed Brown   }
73c4762a1bSJed Brown   PetscFunctionReturn(0);
74c4762a1bSJed Brown }
75c4762a1bSJed Brown 
769371c9d4SSatish Balay static PetscErrorCode RHSJacobianP(TS ts, PetscReal t, Vec X, Mat A, void *ctx) {
77c4762a1bSJed Brown   User               user  = (User)ctx;
78c4762a1bSJed Brown   PetscInt           row[] = {0, 1}, col[] = {0};
79c4762a1bSJed Brown   PetscScalar        J[2][1];
80c4762a1bSJed Brown   const PetscScalar *x;
81c4762a1bSJed Brown 
82c4762a1bSJed Brown   PetscFunctionBeginUser;
83c4762a1bSJed Brown   if (user->combined) col[0] = 2;
849566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
85c4762a1bSJed Brown   J[0][0] = 0;
86c4762a1bSJed Brown   J[1][0] = (1. - x[0] * x[0]) * x[1] - x[0];
879566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A, 2, row, 1, col, &J[0][0], INSERT_VALUES));
889566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
89c4762a1bSJed Brown 
909566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY));
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY));
92c4762a1bSJed Brown   PetscFunctionReturn(0);
93c4762a1bSJed Brown }
94c4762a1bSJed Brown 
95c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
969371c9d4SSatish Balay static PetscErrorCode Monitor(TS ts, PetscInt step, PetscReal t, Vec X, void *ctx) {
97c4762a1bSJed Brown   const PetscScalar *x;
98c4762a1bSJed Brown   PetscReal          tfinal, dt;
99c4762a1bSJed Brown   User               user = (User)ctx;
100c4762a1bSJed Brown   Vec                interpolatedX;
101c4762a1bSJed Brown 
102c4762a1bSJed Brown   PetscFunctionBeginUser;
1039566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts, &dt));
1049566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts, &tfinal));
105c4762a1bSJed Brown 
106c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
1079566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X, &interpolatedX));
1089566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts, user->next_output, interpolatedX));
1099566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX, &x));
1109371c9d4SSatish Balay     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])));
1119566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX, &x));
1129566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
113c4762a1bSJed Brown     user->next_output += 0.1;
114c4762a1bSJed Brown   }
115c4762a1bSJed Brown   PetscFunctionReturn(0);
116c4762a1bSJed Brown }
117c4762a1bSJed Brown 
1189371c9d4SSatish Balay int main(int argc, char **argv) {
119c4762a1bSJed Brown   TS             ts;
120c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
121c4762a1bSJed Brown   PetscScalar   *x_ptr;
122c4762a1bSJed Brown   PetscMPIInt    size;
123c4762a1bSJed Brown   struct _n_User user;
124c4762a1bSJed Brown   PetscInt       rows, cols;
125c4762a1bSJed Brown 
126c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
127c4762a1bSJed Brown      Initialize program
128c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
129327415f7SBarry Smith   PetscFunctionBeginUser;
1309566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
131c4762a1bSJed Brown 
1329566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
1333c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136c4762a1bSJed Brown     Set runtime options
137c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138c4762a1bSJed Brown   user.next_output = 0.0;
139c4762a1bSJed Brown   user.mu          = 1.0e6;
140c4762a1bSJed Brown   user.steps       = 0;
141c4762a1bSJed Brown   user.ftime       = 0.5;
142c4762a1bSJed Brown   user.combined    = PETSC_FALSE;
1439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-monitor", &monitor, NULL));
1449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL, NULL, "-mu", &user.mu, NULL));
1459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL, NULL, "-combined", &user.combined, NULL));
146c4762a1bSJed Brown 
147c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
148c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
149c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
150c4762a1bSJed Brown   rows = 2;
151c4762a1bSJed Brown   cols = user.combined ? 3 : 1;
1529566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jac));
1539566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jac, PETSC_DECIDE, PETSC_DECIDE, 2, 2));
1549566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jac));
1559566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jac));
1569566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jac, &user.x, NULL));
157c4762a1bSJed Brown 
158c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159c4762a1bSJed Brown      Create timestepping solver context
160c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1619566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1629566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSBEULER));
1639566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, IFunction, &user));
1649566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, user.Jac, user.Jac, IJacobian, &user));
1659566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
1669566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, user.ftime));
167*48a46eb9SPierre Jolivet   if (monitor) PetscCall(TSMonitorSet(ts, Monitor, &user, NULL));
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170c4762a1bSJed Brown      Set initial conditions
171c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1729566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.x, &x_ptr));
1739371c9d4SSatish Balay   x_ptr[0] = 2.0;
1749371c9d4SSatish Balay   x_ptr[1] = -0.66666654321;
1759566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.x, &x_ptr));
1769566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 1.0 / 1024.0));
177c4762a1bSJed Brown 
178c4762a1bSJed Brown   /* Set up forward sensitivity */
1799566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Jacp));
1809566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jacp, PETSC_DECIDE, PETSC_DECIDE, rows, cols));
1819566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jacp));
1829566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jacp));
1839566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD, PETSC_DECIDE, PETSC_DECIDE, rows, cols, NULL, &user.sp));
184c4762a1bSJed Brown   if (user.combined) {
1859566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(user.sp));
1869566063dSJacob Faibussowitsch     PetscCall(MatShift(user.sp, 1.0));
187c4762a1bSJed Brown   } else {
1889566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(user.sp));
189c4762a1bSJed Brown   }
1909566063dSJacob Faibussowitsch   PetscCall(TSForwardSetSensitivities(ts, cols, user.sp));
1919566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts, user.Jacp, RHSJacobianP, &user));
192c4762a1bSJed Brown 
193c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194c4762a1bSJed Brown      Set runtime options
195c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1969566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
197c4762a1bSJed Brown 
1989566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, user.x));
1999566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &user.ftime));
2009566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &user.steps));
20163a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "mu %g, steps %" PetscInt_FMT ", ftime %g\n", (double)user.mu, user.steps, (double)user.ftime));
2029566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n ode solution \n"));
2039566063dSJacob Faibussowitsch   PetscCall(VecView(user.x, PETSC_VIEWER_STDOUT_WORLD));
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   if (user.combined) {
2069566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n"));
207c4762a1bSJed Brown   } else {
2089566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n forward sensitivity: d[y(tf) z(tf)]/d[mu]\n"));
209c4762a1bSJed Brown   }
2109566063dSJacob Faibussowitsch   PetscCall(MatView(user.sp, PETSC_VIEWER_STDOUT_WORLD));
211c4762a1bSJed Brown 
212c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
214c4762a1bSJed Brown      are no longer needed.
215c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2169566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jac));
2179566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.sp));
2189566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jacp));
2199566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.x));
2209566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
221c4762a1bSJed Brown 
2229566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
223b122ec5aSJacob Faibussowitsch   return 0;
224c4762a1bSJed Brown }
225c4762a1bSJed Brown 
226c4762a1bSJed Brown /*TEST
227c4762a1bSJed Brown 
228c4762a1bSJed Brown     test:
229c4762a1bSJed Brown       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -combined
230c4762a1bSJed Brown       requires:  !complex !single
231c4762a1bSJed Brown 
232c4762a1bSJed Brown     test:
233c4762a1bSJed Brown       suffix: 2
234c4762a1bSJed Brown       requires: !complex !single
235c4762a1bSJed Brown       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5
236c4762a1bSJed Brown 
237c4762a1bSJed Brown TEST*/
238