Actual source code: ex20adj.c

petsc-3.7.3 2016-08-01
Report Typos and Errors
  1: #define c11 1.0
  2: #define c12 0
  3: #define c21 2.0
  4: #define c22 1.0
  5: static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n";

  7: /*
  8:    Concepts: TS^time-dependent nonlinear problems
  9:    Concepts: TS^van der Pol equation DAE equivalent
 10:    Concepts: TS^adjoint sensitivity analysis
 11:    Processors: 1
 12: */
 13: /* ------------------------------------------------------------------------

 15:    This program solves the van der Pol DAE ODE equivalent
 16:       [ u_1' ] = [      u_2              ]  (2)
 17:       [ u_2' ]   [ mu[(1-u_1^2)u_2-u_1]  ]
 18:    on the domain 0 <= x <= 1, with the boundary conditions
 19:        u_1(0) = 2, u_2(0) = -6.666665432100101e-01,
 20:    and
 21:        mu = 10^6,
 22:    and computes the sensitivities of the final solution w.r.t. initial conditions and parameter \mu with the implicit theta method and its discrete adjoint.

 24:    Notes:
 25:    This code demonstrates the TSAdjoint interface to a DAE system.

 27:    The user provides the implicit right-hand-side function
 28:    [ G(u',u,t) ] = [u' - f(u,t)] = [ u_1'] - [        u_2           ]
 29:                                    [ u_2']   [ mu[(1-u_1^2)u_2-u_1] ]

 31:    and the Jacobian of G (from the PETSc user manual)

 33:               dG   dG
 34:    J(G) = a * -- + --
 35:               du'  du

 37:    and the JacobianP of the explicit right-hand side of (2) f(u,t) ( which is equivalent to -G(0,u,t) ).
 38:    df   [       0         ]
 39:    -- = [                 ]
 40:    dp   [ (1 - u_1^2) u_2 ].

 42:    See ex20.c for more details on the Jacobian.

 44:    Many DAEs can be represented in a general form M u_t = f(u,t).
 45:    Thus both sides of (1) are multiplied by an artificial matrix
 46:    M = [ c11 c12 ]
 47:        [ c21 c22 ]
 48:    to turn (1) into the general form. This operation does not change the solution and it is intended for illustration only.

 50:   ------------------------------------------------------------------------- */
 51: #include <petscts.h>
 52: #include <petsctao.h>

 54: typedef struct _n_User *User;
 55: struct _n_User {
 56:   PetscReal mu;
 57:   PetscReal next_output;

 59:   /* Sensitivity analysis support */
 60:   PetscInt  steps;
 61:   PetscReal ftime;
 62:   Mat       A;                       /* Jacobian matrix */
 63:   Mat       Jacp;                    /* JacobianP matrix */
 64:   Vec       x,lambda[2],mup[2];  /* adjoint variables */
 65: };

 67: /*
 68: *  User-defined routines
 69: */
 72: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 73: {
 74:   PetscErrorCode    ierr;
 75:   User              user = (User)ctx;
 76:   const PetscScalar *x,*xdot;
 77:   PetscScalar       *f;

 80:   VecGetArrayRead(X,&x);
 81:   VecGetArrayRead(Xdot,&xdot);
 82:   VecGetArray(F,&f);
 83:   f[0] = xdot[0] - x[1];
 84:   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
 85:   VecRestoreArrayRead(X,&x);
 86:   VecRestoreArrayRead(Xdot,&xdot);
 87:   VecRestoreArray(F,&f);
 88:   return(0);
 89: }

 93: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 94: {
 95:   PetscErrorCode    ierr;
 96:   User              user     = (User)ctx;
 97:   PetscInt          rowcol[] = {0,1};
 98:   PetscScalar       J[2][2];
 99:   const PetscScalar *x;

102:   VecGetArrayRead(X,&x);

104:   J[0][0] = a;     J[0][1] =  -1.0;
105:   J[1][0] = c21*a + user->mu*(1.0 + 2.0*x[0]*x[1]);   J[1][1] = -c21 + a - user->mu*(1.0-x[0]*x[0]);

107:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
108:   VecRestoreArrayRead(X,&x);

110:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
111:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
112:   if (A != B) {
113:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
114:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
115:   }
116:   return(0);
117: }

121: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
122: {
123:   PetscErrorCode    ierr;
124:   PetscInt          row[] = {0,1},col[]={0};
125:   PetscScalar       J[2][1];
126:   const PetscScalar *x;

129:   VecGetArrayRead(X,&x);

131:   J[0][0] = 0;
132:   J[1][0] = (1.-x[0]*x[0])*x[1]-x[0];
133:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

135:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
136:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
137:   return(0);
138: }

142: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
143: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
144: {
145:   PetscErrorCode    ierr;
146:   const PetscScalar *x;
147:   PetscReal         tfinal, dt;
148:   User              user = (User)ctx;
149:   Vec               interpolatedX;

152:   TSGetTimeStep(ts,&dt);
153:   TSGetDuration(ts,NULL,&tfinal);

155:   while (user->next_output <= t && user->next_output <= tfinal) {
156:     VecDuplicate(X,&interpolatedX);
157:     TSInterpolate(ts,user->next_output,interpolatedX);
158:     VecGetArrayRead(interpolatedX,&x);
159:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
160:                        user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
161:                        (double)PetscRealPart(x[1]));
162:     VecRestoreArrayRead(interpolatedX,&x);
163:     VecDestroy(&interpolatedX);
164:     user->next_output += 0.1;
165:   }
166:   return(0);
167: }

171: int main(int argc,char **argv)
172: {
173:   TS             ts;            /* nonlinear solver */
174:   PetscBool      monitor = PETSC_FALSE;
175:   PetscScalar    *x_ptr,*y_ptr;
176:   PetscMPIInt    size;
177:   struct _n_User user;

180:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181:      Initialize program
182:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183:   PetscInitialize(&argc,&argv,NULL,help);

185:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
186:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

188:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189:     Set runtime options
190:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191:   user.next_output = 0.0;
192:   user.mu          = 1.0e6;
193:   user.steps       = 0;
194:   user.ftime       = 0.5;
195:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
196:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);

198:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199:     Create necessary matrix and vectors, solve same ODE on every process
200:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201:   MatCreate(PETSC_COMM_WORLD,&user.A);
202:   MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
203:   MatSetFromOptions(user.A);
204:   MatSetUp(user.A);
205:   MatCreateVecs(user.A,&user.x,NULL);

207:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
208:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
209:   MatSetFromOptions(user.Jacp);
210:   MatSetUp(user.Jacp);

212:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
213:      Create timestepping solver context
214:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
215:   TSCreate(PETSC_COMM_WORLD,&ts);
216:   TSSetType(ts,TSCN);
217:   TSSetIFunction(ts,NULL,IFunction,&user);
218:   TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);
219:   TSSetDuration(ts,200000,user.ftime);
220:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
221:   if (monitor) {
222:     TSMonitorSet(ts,Monitor,&user,NULL);
223:   }

225:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226:      Set initial conditions
227:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
228:   VecGetArray(user.x,&x_ptr);
229:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
230:   VecRestoreArray(user.x,&x_ptr);
231:   TSSetInitialTimeStep(ts,0.0,.0001);

233:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234:     Save trajectory of solution so that TSAdjointSolve() may be used
235:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236:   TSSetSaveTrajectory(ts);

238:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239:      Set runtime options
240:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
241:   TSSetFromOptions(ts);

243:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
244:      Solve nonlinear system
245:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
246:   TSSolve(ts,user.x);
247:   TSGetSolveTime(ts,&user.ftime);
248:   TSGetTimeStepNumber(ts,&user.steps);

250:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
251:      Adjoint model starts here
252:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
253:   MatCreateVecs(user.A,&user.lambda[0],NULL);
254:   /*   Set initial conditions for the adjoint integration */
255:   VecGetArray(user.lambda[0],&y_ptr);
256:   y_ptr[0] = 1.0; y_ptr[1] = 0.0;
257:   VecRestoreArray(user.lambda[0],&y_ptr);
258:   MatCreateVecs(user.A,&user.lambda[1],NULL);
259:   VecGetArray(user.lambda[1],&y_ptr);
260:   y_ptr[0] = 0.0; y_ptr[1] = 1.0;
261:   VecRestoreArray(user.lambda[1],&y_ptr);

263:   MatCreateVecs(user.Jacp,&user.mup[0],NULL);
264:   VecGetArray(user.mup[0],&x_ptr);
265:   x_ptr[0] = 0.0;
266:   VecRestoreArray(user.mup[0],&x_ptr);
267:   MatCreateVecs(user.Jacp,&user.mup[1],NULL);
268:   VecGetArray(user.mup[1],&x_ptr);
269:   x_ptr[0] = 0.0;
270:   VecRestoreArray(user.mup[1],&x_ptr);

272:   TSSetCostGradients(ts,2,user.lambda,user.mup);

274:   /*   Set RHS JacobianP */
275:   TSAdjointSetRHSJacobian(ts,user.Jacp,RHSJacobianP,&user);

277:   TSAdjointSolve(ts);

279:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[y(tf)]/d[y0]  d[y(tf)]/d[z0]\n");
280:   VecView(user.lambda[0],PETSC_VIEWER_STDOUT_WORLD);
281:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[z(tf)]/d[y0]  d[z(tf)]/d[z0]\n");
282:   VecView(user.lambda[1],PETSC_VIEWER_STDOUT_WORLD);
283:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameters: d[y(tf)]/d[mu]\n");
284:   VecView(user.mup[0],PETSC_VIEWER_STDOUT_WORLD);
285:   PetscPrintf(PETSC_COMM_WORLD,"\n sensivitity wrt parameters: d[z(tf)]/d[mu]\n");
286:   VecView(user.mup[1],PETSC_VIEWER_STDOUT_WORLD);

288:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
289:      Free work space.  All PETSc objects should be destroyed when they
290:      are no longer needed.
291:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
292:   MatDestroy(&user.A);
293:   MatDestroy(&user.Jacp);
294:   VecDestroy(&user.x);
295:   VecDestroy(&user.lambda[0]);
296:   VecDestroy(&user.lambda[1]);
297:   VecDestroy(&user.mup[0]);
298:   VecDestroy(&user.mup[1]);
299:   TSDestroy(&ts);

301:   PetscFinalize();
302:   return(0);
303: }