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