xref: /petsc/src/ts/tutorials/ex23fwdadj.c (revision 82ad9101c7b8cc0cb17c066decc247deb4a8323b)
1*82ad9101SHong Zhang static char help[] = "A toy example for testing forward and adjoint sensitivity analysis of an implicit ODE with a paramerized mass matrice.\n";
2*82ad9101SHong Zhang 
3*82ad9101SHong Zhang /*
4*82ad9101SHong Zhang   This example solves the simple ODE
5*82ad9101SHong Zhang     c x' = b x, x(0) = a,
6*82ad9101SHong Zhang   whose analytical solution is x(T)=a*exp(b/c*T), and calculates the derivative of x(T) w.r.t. c (by default) or w.r.t. b (can be enabled with command line option -der 2).
7*82ad9101SHong Zhang 
8*82ad9101SHong Zhang */
9*82ad9101SHong Zhang 
10*82ad9101SHong Zhang #include <petscts.h>
11*82ad9101SHong Zhang 
12*82ad9101SHong Zhang typedef struct _n_User *User;
13*82ad9101SHong Zhang struct _n_User {
14*82ad9101SHong Zhang   PetscReal a;
15*82ad9101SHong Zhang   PetscReal b;
16*82ad9101SHong Zhang   PetscReal c;
17*82ad9101SHong Zhang   /* Sensitivity analysis support */
18*82ad9101SHong Zhang   PetscInt  steps;
19*82ad9101SHong Zhang   PetscReal ftime;
20*82ad9101SHong Zhang   Mat       Jac;                    /* Jacobian matrix */
21*82ad9101SHong Zhang   Mat       Jacp;                   /* JacobianP matrix */
22*82ad9101SHong Zhang   Vec       x;
23*82ad9101SHong Zhang   Mat       sp;                     /* forward sensitivity variables */
24*82ad9101SHong Zhang   Vec       lambda[1];              /* adjoint sensitivity variables */
25*82ad9101SHong Zhang   Vec       mup[1];                 /* adjoint sensitivity variables */
26*82ad9101SHong Zhang   PetscInt  der;
27*82ad9101SHong Zhang };
28*82ad9101SHong Zhang 
29*82ad9101SHong Zhang static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
30*82ad9101SHong Zhang {
31*82ad9101SHong Zhang   PetscErrorCode    ierr;
32*82ad9101SHong Zhang   User              user = (User)ctx;
33*82ad9101SHong Zhang   const PetscScalar *x,*xdot;
34*82ad9101SHong Zhang   PetscScalar       *f;
35*82ad9101SHong Zhang 
36*82ad9101SHong Zhang   PetscFunctionBeginUser;
37*82ad9101SHong Zhang   ierr = VecGetArrayRead(X,&x);CHKERRQ(ierr);
38*82ad9101SHong Zhang   ierr = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
39*82ad9101SHong Zhang   ierr = VecGetArrayWrite(F,&f);CHKERRQ(ierr);
40*82ad9101SHong Zhang   f[0] = user->c*xdot[0] - user->b*x[0];
41*82ad9101SHong Zhang   ierr = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
42*82ad9101SHong Zhang   ierr = VecRestoreArrayRead(Xdot,&xdot);CHKERRQ(ierr);
43*82ad9101SHong Zhang   ierr = VecRestoreArrayWrite(F,&f);CHKERRQ(ierr);
44*82ad9101SHong Zhang   PetscFunctionReturn(0);
45*82ad9101SHong Zhang }
46*82ad9101SHong Zhang 
47*82ad9101SHong Zhang static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
48*82ad9101SHong Zhang {
49*82ad9101SHong Zhang   PetscErrorCode    ierr;
50*82ad9101SHong Zhang   User              user     = (User)ctx;
51*82ad9101SHong Zhang   PetscInt          rowcol[] = {0};
52*82ad9101SHong Zhang   PetscScalar       J[1][1];
53*82ad9101SHong Zhang   const PetscScalar *x;
54*82ad9101SHong Zhang 
55*82ad9101SHong Zhang   PetscFunctionBeginUser;
56*82ad9101SHong Zhang   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
57*82ad9101SHong Zhang   J[0][0] = user->c*a - user->b*1.0;
58*82ad9101SHong Zhang   ierr    = MatSetValues(B,1,rowcol,1,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
59*82ad9101SHong Zhang   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
60*82ad9101SHong Zhang 
61*82ad9101SHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
62*82ad9101SHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
63*82ad9101SHong Zhang   if (A != B) {
64*82ad9101SHong Zhang     ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
65*82ad9101SHong Zhang     ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
66*82ad9101SHong Zhang   }
67*82ad9101SHong Zhang   PetscFunctionReturn(0);
68*82ad9101SHong Zhang }
69*82ad9101SHong Zhang 
70*82ad9101SHong Zhang static PetscErrorCode IJacobianP(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal shift,Mat A,void *ctx)
71*82ad9101SHong Zhang {
72*82ad9101SHong Zhang   User              user = (User)ctx;
73*82ad9101SHong Zhang   PetscInt          row[] = {0},col[]={0};
74*82ad9101SHong Zhang   PetscScalar       J[1][1];
75*82ad9101SHong Zhang   const PetscScalar *x,*xdot;
76*82ad9101SHong Zhang   PetscReal         dt;
77*82ad9101SHong Zhang   PetscErrorCode    ierr;
78*82ad9101SHong Zhang 
79*82ad9101SHong Zhang   PetscFunctionBeginUser;
80*82ad9101SHong Zhang   ierr    = VecGetArrayRead(X,&x);CHKERRQ(ierr);
81*82ad9101SHong Zhang   ierr    = VecGetArrayRead(Xdot,&xdot);CHKERRQ(ierr);
82*82ad9101SHong Zhang   ierr    = TSGetTimeStep(ts,&dt);CHKERRQ(ierr);
83*82ad9101SHong Zhang   if (user->der == 1) J[0][0] = xdot[0];
84*82ad9101SHong Zhang   if (user->der == 2) J[0][0] = -x[0];
85*82ad9101SHong Zhang   ierr    = MatSetValues(A,1,row,1,col,&J[0][0],INSERT_VALUES);CHKERRQ(ierr);
86*82ad9101SHong Zhang   ierr    = VecRestoreArrayRead(X,&x);CHKERRQ(ierr);
87*82ad9101SHong Zhang 
88*82ad9101SHong Zhang   ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
89*82ad9101SHong Zhang   ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
90*82ad9101SHong Zhang   PetscFunctionReturn(0);
91*82ad9101SHong Zhang }
92*82ad9101SHong Zhang 
93*82ad9101SHong Zhang int main(int argc,char **argv)
94*82ad9101SHong Zhang {
95*82ad9101SHong Zhang   TS             ts;
96*82ad9101SHong Zhang   PetscScalar    *x_ptr;
97*82ad9101SHong Zhang   PetscMPIInt    size;
98*82ad9101SHong Zhang   struct _n_User user;
99*82ad9101SHong Zhang   PetscInt       rows,cols;
100*82ad9101SHong Zhang   PetscErrorCode ierr;
101*82ad9101SHong Zhang 
102*82ad9101SHong Zhang   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
103*82ad9101SHong Zhang 
104*82ad9101SHong Zhang   ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRQ(ierr);
105*82ad9101SHong Zhang   if (size != 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
106*82ad9101SHong Zhang 
107*82ad9101SHong Zhang   user.a           = 2.0;
108*82ad9101SHong Zhang   user.b           = 4.0;
109*82ad9101SHong Zhang   user.c           = 3.0;
110*82ad9101SHong Zhang   user.steps       = 0;
111*82ad9101SHong Zhang   user.ftime       = 1.0;
112*82ad9101SHong Zhang   user.der         = 1;
113*82ad9101SHong Zhang   ierr = PetscOptionsGetInt(NULL,NULL,"-der",&user.der,NULL);CHKERRQ(ierr);
114*82ad9101SHong Zhang 
115*82ad9101SHong Zhang   rows = 1;
116*82ad9101SHong Zhang   cols = 1;
117*82ad9101SHong Zhang   ierr = MatCreate(PETSC_COMM_WORLD,&user.Jac);CHKERRQ(ierr);
118*82ad9101SHong Zhang   ierr = MatSetSizes(user.Jac,PETSC_DECIDE,PETSC_DECIDE,1,1);CHKERRQ(ierr);
119*82ad9101SHong Zhang   ierr = MatSetFromOptions(user.Jac);CHKERRQ(ierr);
120*82ad9101SHong Zhang   ierr = MatSetUp(user.Jac);CHKERRQ(ierr);
121*82ad9101SHong Zhang   ierr = MatCreateVecs(user.Jac,&user.x,NULL);CHKERRQ(ierr);
122*82ad9101SHong Zhang 
123*82ad9101SHong Zhang   ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr);
124*82ad9101SHong Zhang   ierr = TSSetType(ts,TSBEULER);CHKERRQ(ierr);
125*82ad9101SHong Zhang   ierr = TSSetIFunction(ts,NULL,IFunction,&user);CHKERRQ(ierr);
126*82ad9101SHong Zhang   ierr = TSSetIJacobian(ts,user.Jac,user.Jac,IJacobian,&user);CHKERRQ(ierr);
127*82ad9101SHong Zhang   ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr);
128*82ad9101SHong Zhang   ierr = TSSetMaxTime(ts,user.ftime);CHKERRQ(ierr);
129*82ad9101SHong Zhang 
130*82ad9101SHong Zhang   ierr = VecGetArrayWrite(user.x,&x_ptr);CHKERRQ(ierr);
131*82ad9101SHong Zhang   x_ptr[0] = user.a;
132*82ad9101SHong Zhang   ierr = VecRestoreArrayWrite(user.x,&x_ptr);CHKERRQ(ierr);
133*82ad9101SHong Zhang   ierr = TSSetTimeStep(ts,0.001);CHKERRQ(ierr);
134*82ad9101SHong Zhang 
135*82ad9101SHong Zhang   /* Set up forward sensitivity */
136*82ad9101SHong Zhang   ierr = MatCreate(PETSC_COMM_WORLD,&user.Jacp);CHKERRQ(ierr);
137*82ad9101SHong Zhang   ierr = MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,rows,cols);CHKERRQ(ierr);
138*82ad9101SHong Zhang   ierr = MatSetFromOptions(user.Jacp);CHKERRQ(ierr);
139*82ad9101SHong Zhang   ierr = MatSetUp(user.Jacp);CHKERRQ(ierr);
140*82ad9101SHong Zhang   ierr = MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,rows,cols,NULL,&user.sp);CHKERRQ(ierr);
141*82ad9101SHong Zhang   ierr = MatZeroEntries(user.sp);CHKERRQ(ierr);
142*82ad9101SHong Zhang   ierr = TSForwardSetSensitivities(ts,cols,user.sp);CHKERRQ(ierr);
143*82ad9101SHong Zhang   ierr = TSSetIJacobianP(ts,user.Jacp,IJacobianP,&user);CHKERRQ(ierr);
144*82ad9101SHong Zhang 
145*82ad9101SHong Zhang   ierr = TSSetSaveTrajectory(ts);CHKERRQ(ierr);
146*82ad9101SHong Zhang   ierr = TSSetFromOptions(ts);CHKERRQ(ierr);
147*82ad9101SHong Zhang 
148*82ad9101SHong Zhang   ierr = TSSolve(ts,user.x);CHKERRQ(ierr);
149*82ad9101SHong Zhang   ierr = TSGetSolveTime(ts,&user.ftime);CHKERRQ(ierr);
150*82ad9101SHong Zhang   ierr = TSGetStepNumber(ts,&user.steps);CHKERRQ(ierr);
151*82ad9101SHong Zhang   ierr = VecGetArray(user.x,&x_ptr);CHKERRQ(ierr);
152*82ad9101SHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n ode solution %g\n",(double)PetscRealPart(x_ptr[0]));CHKERRQ(ierr);
153*82ad9101SHong Zhang   ierr = VecRestoreArray(user.x,&x_ptr);CHKERRQ(ierr);
154*82ad9101SHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n analytical solution %g\n",(double)user.a*PetscExpReal(user.b/user.c*user.ftime));CHKERRQ(ierr);
155*82ad9101SHong Zhang 
156*82ad9101SHong Zhang   if (user.der == 1) {
157*82ad9101SHong Zhang     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n analytical derivative w.r.t. c %g\n",(double)-user.a*user.ftime*user.b/(user.c*user.c)*PetscExpReal(user.b/user.c*user.ftime));CHKERRQ(ierr);
158*82ad9101SHong Zhang   }
159*82ad9101SHong Zhang   if (user.der == 2) {
160*82ad9101SHong Zhang     ierr = PetscPrintf(PETSC_COMM_WORLD,"\n analytical derivative w.r.t. b %g\n",user.a*user.ftime/user.c*PetscExpReal(user.b/user.c*user.ftime));CHKERRQ(ierr);
161*82ad9101SHong Zhang   }
162*82ad9101SHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n forward sensitivity:\n");CHKERRQ(ierr);
163*82ad9101SHong Zhang   ierr = MatView(user.sp,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
164*82ad9101SHong Zhang 
165*82ad9101SHong Zhang   ierr = MatCreateVecs(user.Jac,&user.lambda[0],NULL);CHKERRQ(ierr);
166*82ad9101SHong Zhang   /* Set initial conditions for the adjoint integration */
167*82ad9101SHong Zhang   ierr = VecGetArrayWrite(user.lambda[0],&x_ptr);CHKERRQ(ierr);
168*82ad9101SHong Zhang   x_ptr[0] = 1.0;
169*82ad9101SHong Zhang   ierr = VecRestoreArrayWrite(user.lambda[0],&x_ptr);CHKERRQ(ierr);
170*82ad9101SHong Zhang   ierr = MatCreateVecs(user.Jacp,&user.mup[0],NULL);CHKERRQ(ierr);
171*82ad9101SHong Zhang   ierr = VecGetArrayWrite(user.mup[0],&x_ptr);CHKERRQ(ierr);
172*82ad9101SHong Zhang   x_ptr[0] = 0.0;
173*82ad9101SHong Zhang   ierr = VecRestoreArrayWrite(user.mup[0],&x_ptr);CHKERRQ(ierr);
174*82ad9101SHong Zhang 
175*82ad9101SHong Zhang   ierr = TSSetCostGradients(ts,1,user.lambda,user.mup);CHKERRQ(ierr);
176*82ad9101SHong Zhang   ierr = TSAdjointSolve(ts);CHKERRQ(ierr);
177*82ad9101SHong Zhang 
178*82ad9101SHong Zhang   ierr = PetscPrintf(PETSC_COMM_WORLD,"\n adjoint sensitivity:\n");CHKERRQ(ierr);
179*82ad9101SHong Zhang   ierr = VecView(user.mup[0],PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
180*82ad9101SHong Zhang 
181*82ad9101SHong Zhang   ierr = MatDestroy(&user.Jac);CHKERRQ(ierr);
182*82ad9101SHong Zhang   ierr = MatDestroy(&user.sp);CHKERRQ(ierr);
183*82ad9101SHong Zhang   ierr = MatDestroy(&user.Jacp);CHKERRQ(ierr);
184*82ad9101SHong Zhang   ierr = VecDestroy(&user.x);CHKERRQ(ierr);
185*82ad9101SHong Zhang   ierr = VecDestroy(&user.lambda[0]);CHKERRQ(ierr);
186*82ad9101SHong Zhang   ierr = VecDestroy(&user.mup[0]);CHKERRQ(ierr);
187*82ad9101SHong Zhang   ierr = TSDestroy(&ts);CHKERRQ(ierr);
188*82ad9101SHong Zhang 
189*82ad9101SHong Zhang   ierr = PetscFinalize();
190*82ad9101SHong Zhang   return(ierr);
191*82ad9101SHong Zhang }
192*82ad9101SHong Zhang 
193*82ad9101SHong Zhang /*TEST
194*82ad9101SHong Zhang 
195*82ad9101SHong Zhang     test:
196*82ad9101SHong Zhang       args: -ts_type beuler
197*82ad9101SHong Zhang 
198*82ad9101SHong Zhang     test:
199*82ad9101SHong Zhang       suffix: 2
200*82ad9101SHong Zhang       args: -ts_type cn
201*82ad9101SHong Zhang       output_file: output/ex23fwdadj_1.out
202*82ad9101SHong Zhang 
203*82ad9101SHong Zhang TEST*/
204