xref: /petsc/src/ts/tutorials/ex16fwd.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
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 #include <petscts.h>
53 #include <petscmat.h>
54 typedef struct _n_User *User;
55 struct _n_User {
56   PetscReal mu;
57   PetscReal next_output;
58   PetscReal tprev;
59 };
60 
61 /*
62    User-defined routines
63 */
64 static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
65 {
66   User              user = (User)ctx;
67   PetscScalar       *f;
68   const PetscScalar *x;
69 
70   PetscFunctionBeginUser;
71   CHKERRQ(VecGetArrayRead(X,&x));
72   CHKERRQ(VecGetArray(F,&f));
73   f[0] = x[1];
74   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
75   CHKERRQ(VecRestoreArrayRead(X,&x));
76   CHKERRQ(VecRestoreArray(F,&f));
77   PetscFunctionReturn(0);
78 }
79 
80 static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
81 {
82   User              user = (User)ctx;
83   PetscReal         mu   = user->mu;
84   PetscInt          rowcol[] = {0,1};
85   PetscScalar       J[2][2];
86   const PetscScalar *x;
87 
88   PetscFunctionBeginUser;
89   CHKERRQ(VecGetArrayRead(X,&x));
90   J[0][0] = 0;
91   J[1][0] = -2.*mu*x[1]*x[0]-1.;
92   J[0][1] = 1.0;
93   J[1][1] = mu*(1.0-x[0]*x[0]);
94   CHKERRQ(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
95   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
96   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
97   if (A != B) {
98     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
99     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
100   }
101   CHKERRQ(VecRestoreArrayRead(X,&x));
102   PetscFunctionReturn(0);
103 }
104 
105 static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
106 {
107   PetscInt          row[] = {0,1},col[]={2};
108   PetscScalar       J[2][1];
109   const PetscScalar *x;
110 
111   PetscFunctionBeginUser;
112   CHKERRQ(VecGetArrayRead(X,&x));
113   J[0][0] = 0;
114   J[1][0] = (1.-x[0]*x[0])*x[1];
115   CHKERRQ(VecRestoreArrayRead(X,&x));
116   CHKERRQ(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
117 
118   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
119   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
120   PetscFunctionReturn(0);
121 }
122 
123 /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
124 static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
125 {
126   const PetscScalar *x;
127   PetscReal         tfinal, dt, tprev;
128   User              user = (User)ctx;
129 
130   PetscFunctionBeginUser;
131   CHKERRQ(TSGetTimeStep(ts,&dt));
132   CHKERRQ(TSGetMaxTime(ts,&tfinal));
133   CHKERRQ(TSGetPrevTime(ts,&tprev));
134   CHKERRQ(VecGetArrayRead(X,&x));
135   CHKERRQ(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])));
136   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev));
137   CHKERRQ(VecRestoreArrayRead(X,&x));
138   PetscFunctionReturn(0);
139 }
140 
141 int main(int argc,char **argv)
142 {
143   TS             ts;            /* nonlinear solver */
144   Vec            x;             /* solution, residual vectors */
145   Mat            A;             /* Jacobian matrix */
146   Mat            Jacp;          /* JacobianP matrix */
147   PetscInt       steps;
148   PetscReal      ftime   =0.5;
149   PetscBool      monitor = PETSC_FALSE;
150   PetscScalar    *x_ptr;
151   PetscMPIInt    size;
152   struct _n_User user;
153   PetscErrorCode ierr;
154   Mat            sp;
155 
156   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
157      Initialize program
158      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
159   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
160   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
161   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
162 
163   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164     Set runtime options
165     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166   user.mu          = 1;
167   user.next_output = 0.0;
168 
169   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
170   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
171 
172   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
173     Create necessary matrix and vectors, solve same ODE on every process
174     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
175   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
176   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
177   CHKERRQ(MatSetFromOptions(A));
178   CHKERRQ(MatSetUp(A));
179   CHKERRQ(MatCreateVecs(A,&x,NULL));
180 
181   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Jacp));
182   CHKERRQ(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,3));
183   CHKERRQ(MatSetFromOptions(Jacp));
184   CHKERRQ(MatSetUp(Jacp));
185 
186   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,3,NULL,&sp));
187   CHKERRQ(MatZeroEntries(sp));
188   CHKERRQ(MatShift(sp,1.0));
189 
190   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191      Create timestepping solver context
192      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
194   CHKERRQ(TSSetType(ts,TSRK));
195   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&user));
196   /*   Set RHS Jacobian for the adjoint integration */
197   CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,&user));
198   CHKERRQ(TSSetMaxTime(ts,ftime));
199   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
200   if (monitor) {
201     CHKERRQ(TSMonitorSet(ts,Monitor,&user,NULL));
202   }
203   CHKERRQ(TSForwardSetSensitivities(ts,3,sp));
204   CHKERRQ(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user));
205 
206   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207      Set initial conditions
208    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209   CHKERRQ(VecGetArray(x,&x_ptr));
210 
211   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
212   CHKERRQ(VecRestoreArray(x,&x_ptr));
213   CHKERRQ(TSSetTimeStep(ts,.001));
214 
215   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216      Set runtime options
217    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
218   CHKERRQ(TSSetFromOptions(ts));
219 
220   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221      Solve nonlinear system
222      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
223   CHKERRQ(TSSolve(ts,x));
224   CHKERRQ(TSGetSolveTime(ts,&ftime));
225   CHKERRQ(TSGetStepNumber(ts,&steps));
226   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime));
227   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
228 
229   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n"));
230   CHKERRQ(MatView(sp,PETSC_VIEWER_STDOUT_WORLD));
231 
232   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233      Free work space.  All PETSc objects should be destroyed when they
234      are no longer needed.
235    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236   CHKERRQ(MatDestroy(&A));
237   CHKERRQ(MatDestroy(&Jacp));
238   CHKERRQ(VecDestroy(&x));
239   CHKERRQ(MatDestroy(&sp));
240   CHKERRQ(TSDestroy(&ts));
241   ierr = PetscFinalize();
242   return ierr;
243 }
244 
245 /*TEST
246 
247     test:
248       args: -monitor 0 -ts_adapt_type none
249 
250 TEST*/
251