xref: /petsc/src/ts/tutorials/ex20fwd.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
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 */
35c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
36c4762a1bSJed Brown {
37c4762a1bSJed Brown   User              user = (User)ctx;
38c4762a1bSJed Brown   const PetscScalar *x,*xdot;
39c4762a1bSJed Brown   PetscScalar       *f;
40c4762a1bSJed Brown 
41c4762a1bSJed Brown   PetscFunctionBeginUser;
429566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
439566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(Xdot,&xdot));
449566063dSJacob Faibussowitsch   PetscCall(VecGetArray(F,&f));
45c4762a1bSJed Brown   f[0] = xdot[0] - x[1];
46c4762a1bSJed Brown   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
479566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(Xdot,&xdot));
499566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(F,&f));
50c4762a1bSJed Brown   PetscFunctionReturn(0);
51c4762a1bSJed Brown }
52c4762a1bSJed Brown 
53c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
54c4762a1bSJed Brown {
55c4762a1bSJed Brown   User              user     = (User)ctx;
56c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
57c4762a1bSJed Brown   PetscScalar       J[2][2];
58c4762a1bSJed Brown   const PetscScalar *x;
59c4762a1bSJed Brown 
60c4762a1bSJed Brown   PetscFunctionBeginUser;
619566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
62c4762a1bSJed Brown   J[0][0] = a;     J[0][1] =  -1.0;
63c4762a1bSJed Brown   J[1][0] = c21*a + user->mu*(1.0 + 2.0*x[0]*x[1]);   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 
76c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
77c4762a1bSJed Brown {
78c4762a1bSJed Brown   User              user = (User)ctx;
79c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[]={0};
80c4762a1bSJed Brown   PetscScalar       J[2][1];
81c4762a1bSJed Brown   const PetscScalar *x;
82c4762a1bSJed Brown 
83c4762a1bSJed Brown   PetscFunctionBeginUser;
84c4762a1bSJed Brown   if (user->combined) col[0] = 2;
859566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X,&x));
86c4762a1bSJed Brown   J[0][0] = 0;
87c4762a1bSJed Brown   J[1][0] = (1.-x[0]*x[0])*x[1]-x[0];
889566063dSJacob Faibussowitsch   PetscCall(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
899566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X,&x));
90c4762a1bSJed Brown 
919566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
929566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
93c4762a1bSJed Brown   PetscFunctionReturn(0);
94c4762a1bSJed Brown }
95c4762a1bSJed Brown 
96c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
97c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
98c4762a1bSJed Brown {
99c4762a1bSJed Brown   const PetscScalar *x;
100c4762a1bSJed Brown   PetscReal         tfinal, dt;
101c4762a1bSJed Brown   User              user = (User)ctx;
102c4762a1bSJed Brown   Vec               interpolatedX;
103c4762a1bSJed Brown 
104c4762a1bSJed Brown   PetscFunctionBeginUser;
1059566063dSJacob Faibussowitsch   PetscCall(TSGetTimeStep(ts,&dt));
1069566063dSJacob Faibussowitsch   PetscCall(TSGetMaxTime(ts,&tfinal));
107c4762a1bSJed Brown 
108c4762a1bSJed Brown   while (user->next_output <= t && user->next_output <= tfinal) {
1099566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X,&interpolatedX));
1109566063dSJacob Faibussowitsch     PetscCall(TSInterpolate(ts,user->next_output,interpolatedX));
1119566063dSJacob Faibussowitsch     PetscCall(VecGetArrayRead(interpolatedX,&x));
11263a3b9bcSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %" PetscInt_FMT " TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
11363a3b9bcSJacob Faibussowitsch                           (double)user->next_output,step,(double)t,(double)dt,
11463a3b9bcSJacob Faibussowitsch                           (double)PetscRealPart(x[0]),(double)PetscRealPart(x[1])));
1159566063dSJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(interpolatedX,&x));
1169566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&interpolatedX));
117c4762a1bSJed Brown     user->next_output += 0.1;
118c4762a1bSJed Brown   }
119c4762a1bSJed Brown   PetscFunctionReturn(0);
120c4762a1bSJed Brown }
121c4762a1bSJed Brown 
122c4762a1bSJed Brown int main(int argc,char **argv)
123c4762a1bSJed Brown {
124c4762a1bSJed Brown   TS             ts;
125c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
126c4762a1bSJed Brown   PetscScalar    *x_ptr;
127c4762a1bSJed Brown   PetscMPIInt    size;
128c4762a1bSJed Brown   struct _n_User user;
129c4762a1bSJed Brown   PetscInt       rows,cols;
130c4762a1bSJed Brown 
131c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132c4762a1bSJed Brown      Initialize program
133c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134*327415f7SBarry Smith   PetscFunctionBeginUser;
1359566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,NULL,help));
136c4762a1bSJed Brown 
1379566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1383c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
141c4762a1bSJed Brown     Set runtime options
142c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
143c4762a1bSJed Brown   user.next_output = 0.0;
144c4762a1bSJed Brown   user.mu          = 1.0e6;
145c4762a1bSJed Brown   user.steps       = 0;
146c4762a1bSJed Brown   user.ftime       = 0.5;
147c4762a1bSJed Brown   user.combined    = PETSC_FALSE;
1489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
1499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
1509566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetBool(NULL,NULL,"-combined",&user.combined,NULL));
151c4762a1bSJed Brown 
152c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
153c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
154c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
155c4762a1bSJed Brown   rows = 2;
156c4762a1bSJed Brown   cols = user.combined ? 3 : 1;
1579566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jac));
1589566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,2,2));
1599566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jac));
1609566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jac));
1619566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(user.Jac,&user.x,NULL));
162c4762a1bSJed Brown 
163c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164c4762a1bSJed Brown      Create timestepping solver context
165c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1669566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD,&ts));
1679566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts,TSBEULER));
1689566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts,NULL,IFunction,&user));
1699566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user));
1709566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
1719566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts,user.ftime));
172c4762a1bSJed Brown   if (monitor) {
1739566063dSJacob Faibussowitsch     PetscCall(TSMonitorSet(ts,Monitor,&user,NULL));
174c4762a1bSJed Brown   }
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177c4762a1bSJed Brown      Set initial conditions
178c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1799566063dSJacob Faibussowitsch   PetscCall(VecGetArray(user.x,&x_ptr));
180c4762a1bSJed Brown   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
1819566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(user.x,&x_ptr));
1829566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts,1.0/1024.0));
183c4762a1bSJed Brown 
184c4762a1bSJed Brown   /* Set up forward sensitivity */
1859566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Jacp));
1869566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols));
1879566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(user.Jacp));
1889566063dSJacob Faibussowitsch   PetscCall(MatSetUp(user.Jacp));
1899566063dSJacob Faibussowitsch   PetscCall(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp));
190c4762a1bSJed Brown   if (user.combined) {
1919566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(user.sp));
1929566063dSJacob Faibussowitsch     PetscCall(MatShift(user.sp,1.0));
193c4762a1bSJed Brown   } else {
1949566063dSJacob Faibussowitsch     PetscCall(MatZeroEntries(user.sp));
195c4762a1bSJed Brown   }
1969566063dSJacob Faibussowitsch   PetscCall(TSForwardSetSensitivities(ts,cols,user.sp));
1979566063dSJacob Faibussowitsch   PetscCall(TSSetRHSJacobianP(ts,user.Jacp,RHSJacobianP,&user));
198c4762a1bSJed Brown 
199c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200c4762a1bSJed Brown      Set runtime options
201c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2029566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
203c4762a1bSJed Brown 
2049566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts,user.x));
2059566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts,&user.ftime));
2069566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts,&user.steps));
20763a3b9bcSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %" PetscInt_FMT ", ftime %g\n",(double)user.mu,user.steps,(double)user.ftime));
2089566063dSJacob Faibussowitsch   PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n ode solution \n"));
2099566063dSJacob Faibussowitsch   PetscCall(VecView(user.x,PETSC_VIEWER_STDOUT_WORLD));
210c4762a1bSJed Brown 
211c4762a1bSJed Brown   if (user.combined) {
2129566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n"));
213c4762a1bSJed Brown   } else {
2149566063dSJacob Faibussowitsch     PetscCall(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[mu]\n"));
215c4762a1bSJed Brown   }
2169566063dSJacob Faibussowitsch   PetscCall(MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD));
217c4762a1bSJed Brown 
218c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
220c4762a1bSJed Brown      are no longer needed.
221c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2229566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jac));
2239566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.sp));
2249566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&user.Jacp));
2259566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&user.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       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -combined
236c4762a1bSJed Brown       requires:  !complex !single
237c4762a1bSJed Brown 
238c4762a1bSJed Brown     test:
239c4762a1bSJed Brown       suffix: 2
240c4762a1bSJed Brown       requires: !complex !single
241c4762a1bSJed Brown       args: -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5
242c4762a1bSJed Brown 
243c4762a1bSJed Brown TEST*/
244