xref: /petsc/src/ts/tutorials/ex16fwd.c (revision b122ec5aa1bd4469eb4e0673542fb7de3f411254)
1c4762a1bSJed Brown static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n\
2c4762a1bSJed Brown Input parameters include:\n\
3c4762a1bSJed Brown       -mu : stiffness parameter\n\n";
4c4762a1bSJed Brown 
5c4762a1bSJed Brown /*
6c4762a1bSJed Brown    Concepts: TS^time-dependent nonlinear problems
7c4762a1bSJed Brown    Concepts: TS^van der Pol equation
8c4762a1bSJed Brown    Concepts: TS^adjoint sensitivity analysis
9c4762a1bSJed Brown    Processors: 1
10c4762a1bSJed Brown */
11c4762a1bSJed Brown /* ------------------------------------------------------------------------
12c4762a1bSJed Brown 
13c4762a1bSJed Brown    This program solves the van der Pol equation
14c4762a1bSJed Brown        y'' - \mu (1-y^2)*y' + y = 0        (1)
15c4762a1bSJed Brown    on the domain 0 <= x <= 1, with the boundary conditions
16c4762a1bSJed Brown        y(0) = 2, y'(0) = 0,
17c4762a1bSJed 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.
18c4762a1bSJed Brown 
19c4762a1bSJed Brown    Notes:
205ab1ac2bSHong Zhang    This code demonstrates the TSForward interface to a system of ordinary differential equations (ODEs) in the form of u_t = f(u,t).
21c4762a1bSJed Brown 
22c4762a1bSJed Brown    (1) can be turned into a system of first order ODEs
23c4762a1bSJed Brown    [ y' ] = [          z          ]
24c4762a1bSJed Brown    [ z' ]   [ \mu (1 - y^2) z - y ]
25c4762a1bSJed Brown 
26c4762a1bSJed Brown    which then we can write as a vector equation
27c4762a1bSJed Brown 
28c4762a1bSJed Brown    [ u_1' ] = [             u_2           ]  (2)
29c4762a1bSJed Brown    [ u_2' ]   [ \mu (1 - u_1^2) u_2 - u_1 ]
30c4762a1bSJed Brown 
31c4762a1bSJed Brown    which is now in the form of u_t = F(u,t).
32c4762a1bSJed Brown 
33c4762a1bSJed Brown    The user provides the right-hand-side function
34c4762a1bSJed Brown 
355ab1ac2bSHong Zhang    [ f(u,t) ] = [ u_2                       ]
36c4762a1bSJed Brown                 [ \mu (1 - u_1^2) u_2 - u_1 ]
37c4762a1bSJed Brown 
38c4762a1bSJed Brown    the Jacobian function
39c4762a1bSJed Brown 
405ab1ac2bSHong Zhang    df   [       0           ;         1        ]
41c4762a1bSJed Brown    -- = [                                      ]
42c4762a1bSJed Brown    du   [ -2 \mu u_1*u_2 - 1;  \mu (1 - u_1^2) ]
43c4762a1bSJed Brown 
44c4762a1bSJed Brown    and the JacobainP (the Jacobian w.r.t. parameter) function
45c4762a1bSJed Brown 
465ab1ac2bSHong Zhang    df      [  0;   0;     0             ]
47c4762a1bSJed Brown    ---   = [                            ]
48c4762a1bSJed Brown    d\mu    [  0;   0;  (1 - u_1^2) u_2  ]
49c4762a1bSJed Brown 
50c4762a1bSJed Brown   ------------------------------------------------------------------------- */
51c4762a1bSJed Brown 
52c4762a1bSJed Brown #include <petscts.h>
53c4762a1bSJed Brown #include <petscmat.h>
54c4762a1bSJed Brown typedef struct _n_User *User;
55c4762a1bSJed Brown struct _n_User {
56c4762a1bSJed Brown   PetscReal mu;
57c4762a1bSJed Brown   PetscReal next_output;
58c4762a1bSJed Brown   PetscReal tprev;
59c4762a1bSJed Brown };
60c4762a1bSJed Brown 
61c4762a1bSJed Brown /*
620e3d61c9SBarry Smith    User-defined routines
63c4762a1bSJed Brown */
64c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
65c4762a1bSJed Brown {
66c4762a1bSJed Brown   User              user = (User)ctx;
67c4762a1bSJed Brown   PetscScalar       *f;
68c4762a1bSJed Brown   const PetscScalar *x;
69c4762a1bSJed Brown 
70c4762a1bSJed Brown   PetscFunctionBeginUser;
715f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
725f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(F,&f));
73c4762a1bSJed Brown   f[0] = x[1];
74c4762a1bSJed Brown   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
755f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
765f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(F,&f));
77c4762a1bSJed Brown   PetscFunctionReturn(0);
78c4762a1bSJed Brown }
79c4762a1bSJed Brown 
80c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
81c4762a1bSJed Brown {
82c4762a1bSJed Brown   User              user = (User)ctx;
83c4762a1bSJed Brown   PetscReal         mu   = user->mu;
84c4762a1bSJed Brown   PetscInt          rowcol[] = {0,1};
85c4762a1bSJed Brown   PetscScalar       J[2][2];
86c4762a1bSJed Brown   const PetscScalar *x;
87c4762a1bSJed Brown 
88c4762a1bSJed Brown   PetscFunctionBeginUser;
895f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
90c4762a1bSJed Brown   J[0][0] = 0;
91c4762a1bSJed Brown   J[1][0] = -2.*mu*x[1]*x[0]-1.;
92c4762a1bSJed Brown   J[0][1] = 1.0;
93c4762a1bSJed Brown   J[1][1] = mu*(1.0-x[0]*x[0]);
945f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES));
955f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
965f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
97c4762a1bSJed Brown   if (A != B) {
985f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY));
995f80ce2aSJacob Faibussowitsch     CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY));
100c4762a1bSJed Brown   }
1015f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
102c4762a1bSJed Brown   PetscFunctionReturn(0);
103c4762a1bSJed Brown }
104c4762a1bSJed Brown 
105c4762a1bSJed Brown static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
106c4762a1bSJed Brown {
107c4762a1bSJed Brown   PetscInt          row[] = {0,1},col[]={2};
108c4762a1bSJed Brown   PetscScalar       J[2][1];
109c4762a1bSJed Brown   const PetscScalar *x;
110c4762a1bSJed Brown 
111c4762a1bSJed Brown   PetscFunctionBeginUser;
1125f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
113c4762a1bSJed Brown   J[0][0] = 0;
114c4762a1bSJed Brown   J[1][0] = (1.-x[0]*x[0])*x[1];
1155f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
1165f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES));
117c4762a1bSJed Brown 
1185f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY));
1195f80ce2aSJacob Faibussowitsch   CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY));
120c4762a1bSJed Brown   PetscFunctionReturn(0);
121c4762a1bSJed Brown }
122c4762a1bSJed Brown 
123c4762a1bSJed Brown /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
124c4762a1bSJed Brown static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
125c4762a1bSJed Brown {
126c4762a1bSJed Brown   const PetscScalar *x;
127c4762a1bSJed Brown   PetscReal         tfinal, dt, tprev;
128c4762a1bSJed Brown   User              user = (User)ctx;
129c4762a1bSJed Brown 
130c4762a1bSJed Brown   PetscFunctionBeginUser;
1315f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetTimeStep(ts,&dt));
1325f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetMaxTime(ts,&tfinal));
1335f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetPrevTime(ts,&tprev));
1345f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArrayRead(X,&x));
1355f80ce2aSJacob Faibussowitsch   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])));
1365f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev));
1375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArrayRead(X,&x));
138c4762a1bSJed Brown   PetscFunctionReturn(0);
139c4762a1bSJed Brown }
140c4762a1bSJed Brown 
141c4762a1bSJed Brown int main(int argc,char **argv)
142c4762a1bSJed Brown {
143c4762a1bSJed Brown   TS             ts;            /* nonlinear solver */
144c4762a1bSJed Brown   Vec            x;             /* solution, residual vectors */
145c4762a1bSJed Brown   Mat            A;             /* Jacobian matrix */
146c4762a1bSJed Brown   Mat            Jacp;          /* JacobianP matrix */
147c4762a1bSJed Brown   PetscInt       steps;
148c4762a1bSJed Brown   PetscReal      ftime   =0.5;
149c4762a1bSJed Brown   PetscBool      monitor = PETSC_FALSE;
150c4762a1bSJed Brown   PetscScalar    *x_ptr;
151c4762a1bSJed Brown   PetscMPIInt    size;
152c4762a1bSJed Brown   struct _n_User user;
153c4762a1bSJed Brown   Mat            sp;
154c4762a1bSJed Brown 
155c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
156c4762a1bSJed Brown      Initialize program
157c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
158*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscInitialize(&argc,&argv,NULL,help));
1595f80ce2aSJacob Faibussowitsch   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
1603c633725SBarry Smith   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
161c4762a1bSJed Brown 
162c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163c4762a1bSJed Brown     Set runtime options
164c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165c4762a1bSJed Brown   user.mu          = 1;
166c4762a1bSJed Brown   user.next_output = 0.0;
167c4762a1bSJed Brown 
1685f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL));
1695f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL));
170c4762a1bSJed Brown 
171c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
172c4762a1bSJed Brown     Create necessary matrix and vectors, solve same ODE on every process
173c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1745f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A));
1755f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2));
1765f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(A));
1775f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(A));
1785f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateVecs(A,&x,NULL));
179c4762a1bSJed Brown 
1805f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&Jacp));
1815f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,3));
1825f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetFromOptions(Jacp));
1835f80ce2aSJacob Faibussowitsch   CHKERRQ(MatSetUp(Jacp));
184c4762a1bSJed Brown 
1855f80ce2aSJacob Faibussowitsch   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,3,NULL,&sp));
1865f80ce2aSJacob Faibussowitsch   CHKERRQ(MatZeroEntries(sp));
1875f80ce2aSJacob Faibussowitsch   CHKERRQ(MatShift(sp,1.0));
188c4762a1bSJed Brown 
189c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190c4762a1bSJed Brown      Create timestepping solver context
191c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1925f80ce2aSJacob Faibussowitsch   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts));
1935f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetType(ts,TSRK));
1945f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSFunction(ts,NULL,RHSFunction,&user));
195c4762a1bSJed Brown   /*   Set RHS Jacobian for the adjoint integration */
1965f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobian(ts,A,A,RHSJacobian,&user));
1975f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetMaxTime(ts,ftime));
1985f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP));
199c4762a1bSJed Brown   if (monitor) {
2005f80ce2aSJacob Faibussowitsch     CHKERRQ(TSMonitorSet(ts,Monitor,&user,NULL));
201c4762a1bSJed Brown   }
2025f80ce2aSJacob Faibussowitsch   CHKERRQ(TSForwardSetSensitivities(ts,3,sp));
2035f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetRHSJacobianP(ts,Jacp,RHSJacobianP,&user));
204c4762a1bSJed Brown 
205c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
206c4762a1bSJed Brown      Set initial conditions
207c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2085f80ce2aSJacob Faibussowitsch   CHKERRQ(VecGetArray(x,&x_ptr));
209c4762a1bSJed Brown 
210c4762a1bSJed Brown   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
2115f80ce2aSJacob Faibussowitsch   CHKERRQ(VecRestoreArray(x,&x_ptr));
2125f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetTimeStep(ts,.001));
213c4762a1bSJed Brown 
214c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215c4762a1bSJed Brown      Set runtime options
216c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2175f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSetFromOptions(ts));
218c4762a1bSJed Brown 
219c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220c4762a1bSJed Brown      Solve nonlinear system
221c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2225f80ce2aSJacob Faibussowitsch   CHKERRQ(TSSolve(ts,x));
2235f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetSolveTime(ts,&ftime));
2245f80ce2aSJacob Faibussowitsch   CHKERRQ(TSGetStepNumber(ts,&steps));
2255f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime));
2265f80ce2aSJacob Faibussowitsch   CHKERRQ(VecView(x,PETSC_VIEWER_STDOUT_WORLD));
227c4762a1bSJed Brown 
2285f80ce2aSJacob Faibussowitsch   CHKERRQ(PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity: d[y(tf) z(tf)]/d[y0 z0 mu]\n"));
2295f80ce2aSJacob Faibussowitsch   CHKERRQ(MatView(sp,PETSC_VIEWER_STDOUT_WORLD));
230c4762a1bSJed Brown 
231c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
232c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they
233c4762a1bSJed Brown      are no longer needed.
234c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2355f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&A));
2365f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&Jacp));
2375f80ce2aSJacob Faibussowitsch   CHKERRQ(VecDestroy(&x));
2385f80ce2aSJacob Faibussowitsch   CHKERRQ(MatDestroy(&sp));
2395f80ce2aSJacob Faibussowitsch   CHKERRQ(TSDestroy(&ts));
240*b122ec5aSJacob Faibussowitsch   CHKERRQ(PetscFinalize());
241*b122ec5aSJacob Faibussowitsch   return 0;
242c4762a1bSJed Brown }
243c4762a1bSJed Brown 
244c4762a1bSJed Brown /*TEST
245c4762a1bSJed Brown 
246c4762a1bSJed Brown     test:
247c4762a1bSJed Brown       args: -monitor 0 -ts_adapt_type none
248c4762a1bSJed Brown 
249c4762a1bSJed Brown TEST*/
250