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