Actual source code: ex20adj.c

petsc-3.6.1 2015-08-06
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[] = "Solves the van der Pol equation.\n\
  6: Input parameters include:\n";

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

 15:    This program solves the van der Pol DAE ODE equivalent
 16:        y' = z                 (1)
 17:        z' = mu[(1-y^2)z-y]
 18:    on the domain 0 <= x <= 1, with the boundary conditions
 19:        y(0) = 2, y'(0) = -6.666665432100101e-01,
 20:    and
 21:        mu = 10^6.
 22:    This is a nonlinear equation.

 24:    Notes:
 25:    This code demonstrates the TS solver interface to a variant of
 26:    linear problems, u_t = f(u,t), namely turning (1) into a system of
 27:    first order differential equations,

 29:    [ y' ] = [          z          ]
 30:    [ z' ]   [     mu[(1-y^2)z-y]  ]

 32:    which then we can write as a vector equation

 34:    [ u_1' ] = [      u_2              ]  (2)
 35:    [ u_2' ]   [ mu[(1-u_1^2)u_2-u_1]  ]

 37:    which is now in the desired form of u_t = f(u,t). One way that we
 38:    can split f(u,t) in (2) is to split by component,

 40:    [ u_1' ] = [  u_2 ] + [       0              ]
 41:    [ u_2' ]   [  0   ]   [ mu[(1-u_1^2)u_2-u_1] ]

 43:    where

 45:    [ F(u,t) ] = [  u_2 ]
 46:                 [  0   ]

 48:    and

 50:    [ G(u',u,t) ] = [ u_1' ] - [            0         ]
 51:                    [ u_2' ]   [ mu[(1-u_1^2)u_2-u_1] ]

 53:    Using the definition of the Jacobian of G (from the PETSc user manual),
 54:    in the equation G(u',u,t) = F(u,t),

 56:               dG   dG
 57:    J(G) = a * -- + --
 58:               du'  du

 60:    where d is the partial derivative. In this example,

 62:    dG   [ 1 ; 0 ]
 63:    -- = [       ]
 64:    du'  [ 0 ; 1 ]

 66:    dG   [ 0                       ;         0         ]
 67:    -- = [                                             ]
 68:    du   [ mu*(1.0 + 2.0*u_1*u_2) ; -mu*(1-u_1*u_1)    ]

 70:    Hence,

 72:           [      a                 ;         0          ]
 73:    J(G) = [                                             ]
 74:           [ mu*(1.0 + 2.0*u_1*u_2) ; a - mu*(1-u_1*u_1) ]

 76:   ------------------------------------------------------------------------- */
 77: #include <petscts.h>
 78: #include <petsctao.h>

 80: typedef struct _n_User *User;
 81: struct _n_User {
 82:   PetscReal mu;
 83:   PetscReal next_output;
 84: 
 85:   /* Sensitivity analysis support */
 86:   PetscInt  steps;
 87:   PetscReal ftime;
 88:   Mat       A;                       /* Jacobian matrix */
 89:   Mat       Jacp;                    /* JacobianP matrix */
 90:   Vec       x,lambda[2],mup[2];  /* adjoint variables */
 91: };

 93: /*
 94: *  User-defined routines
 95: */
 98: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 99: {
100:   PetscErrorCode    ierr;
101:   User              user = (User)ctx;
102:   const PetscScalar *x,*xdot;
103:   PetscScalar       *f;

106:   VecGetArrayRead(X,&x);
107:   VecGetArrayRead(Xdot,&xdot);
108:   VecGetArray(F,&f);
109:   f[0] = xdot[0] - x[1];
110:   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
111:   VecRestoreArrayRead(X,&x);
112:   VecRestoreArrayRead(Xdot,&xdot);
113:   VecRestoreArray(F,&f);
114:   return(0);
115: }

119: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
120: {
121:   PetscErrorCode    ierr;
122:   User              user     = (User)ctx;
123:   PetscInt          rowcol[] = {0,1};
124:   PetscScalar       J[2][2];
125:   const PetscScalar *x;

128:   VecGetArrayRead(X,&x);

130:   J[0][0] = a;     J[0][1] =  -1.0;
131:   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]);
132: 
133:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
134:   VecRestoreArrayRead(X,&x);

136:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
137:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
138:   if (A != B) {
139:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
140:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
141:   }
142:   return(0);
143: }

147: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
148: {
149:   PetscErrorCode    ierr;
150:   PetscInt          row[] = {0,1},col[]={0};
151:   PetscScalar       J[2][1];
152:   const PetscScalar *x;

155:   VecGetArrayRead(X,&x);

157:   J[0][0] = 0;
158:   J[1][0] = (1.-x[0]*x[0])*x[1]-x[0];
159:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

161:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
162:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
163:   return(0);
164: }

168: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
169: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
170: {
171:   PetscErrorCode    ierr;
172:   const PetscScalar *x;
173:   PetscReal         tfinal, dt;
174:   User              user = (User)ctx;
175:   Vec               interpolatedX;

178:   TSGetTimeStep(ts,&dt);
179:   TSGetDuration(ts,NULL,&tfinal);

181:   while (user->next_output <= t && user->next_output <= tfinal) {
182:     VecDuplicate(X,&interpolatedX);
183:     TSInterpolate(ts,user->next_output,interpolatedX);
184:     VecGetArrayRead(interpolatedX,&x);
185:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
186:                        user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
187:                        (double)PetscRealPart(x[1]));
188:     VecRestoreArrayRead(interpolatedX,&x);
189:     VecDestroy(&interpolatedX);
190:     user->next_output += 0.1;
191:   }
192:   return(0);
193: }

197: int main(int argc,char **argv)
198: {
199:   TS             ts;            /* nonlinear solver */
200:   PetscBool      monitor = PETSC_FALSE;
201:   PetscScalar    *x_ptr,*y_ptr;
202:   PetscMPIInt    size;
203:   struct _n_User user;

206:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207:      Initialize program
208:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209:   PetscInitialize(&argc,&argv,NULL,help);
210: 
211:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
212:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

214:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215:     Set runtime options
216:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217:   user.next_output = 0.0;
218:   user.mu          = 1.0e6;
219:   user.steps       = 0;
220:   user.ftime       = 0.5;
221:   PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);
222:   PetscOptionsGetReal(NULL,"-mu",&user.mu,NULL);

224:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
225:     Create necessary matrix and vectors, solve same ODE on every process
226:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
227:   MatCreate(PETSC_COMM_WORLD,&user.A);
228:   MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
229:   MatSetFromOptions(user.A);
230:   MatSetUp(user.A);
231:   MatCreateVecs(user.A,&user.x,NULL);

233:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
234:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
235:   MatSetFromOptions(user.Jacp);
236:   MatSetUp(user.Jacp);

238:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239:      Create timestepping solver context
240:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
241:   TSCreate(PETSC_COMM_WORLD,&ts);
242:   TSSetType(ts,TSCN);
243:   TSSetIFunction(ts,NULL,IFunction,&user);
244:   TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);
245:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
246:   TSSetDuration(ts,200000,user.ftime);
247:   if (monitor) {
248:     TSMonitorSet(ts,Monitor,&user,NULL);
249:   }

251:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
252:      Set initial conditions
253:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
254:   VecGetArray(user.x,&x_ptr);
255:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
256:   VecRestoreArray(user.x,&x_ptr);
257:   TSSetInitialTimeStep(ts,0.0,.0001);

259:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
260:     Save trajectory of solution so that TSAdjointSolve() may be used
261:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
262:   TSSetSaveTrajectory(ts);

264:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265:      Set runtime options
266:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
267:   TSSetFromOptions(ts);

269:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
270:      Solve nonlinear system
271:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
272:   TSSolve(ts,user.x);
273:   TSGetSolveTime(ts,&user.ftime);
274:   TSGetTimeStepNumber(ts,&user.steps);
275:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);
276:   PetscPrintf(PETSC_COMM_WORLD,"\n ode solution \n");
277:   VecView(user.x,PETSC_VIEWER_STDOUT_WORLD);

279:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
280:      Adjoint model starts here
281:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
282:   MatCreateVecs(user.A,&user.lambda[0],NULL);
283:   /*   Set initial conditions for the adjoint integration */
284:   VecGetArray(user.lambda[0],&y_ptr);
285:   y_ptr[0] = 1.0; y_ptr[1] = 0.0;
286:   VecRestoreArray(user.lambda[0],&y_ptr);
287:   MatCreateVecs(user.A,&user.lambda[1],NULL);
288:   VecGetArray(user.lambda[1],&y_ptr);
289:   y_ptr[0] = 0.0; y_ptr[1] = 1.0;
290:   VecRestoreArray(user.lambda[1],&y_ptr);

292:   MatCreateVecs(user.Jacp,&user.mup[0],NULL);
293:   VecGetArray(user.mup[0],&x_ptr);
294:   x_ptr[0] = 0.0;
295:   VecRestoreArray(user.mup[0],&x_ptr);
296:   MatCreateVecs(user.Jacp,&user.mup[1],NULL);
297:   VecGetArray(user.mup[1],&x_ptr);
298:   x_ptr[0] = 0.0;
299:   VecRestoreArray(user.mup[1],&x_ptr);

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

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

306:   TSAdjointSolve(ts);

308:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[y(tf)]/d[y0]  d[y(tf)]/d[z0]\n");
309:   VecView(user.lambda[0],PETSC_VIEWER_STDOUT_WORLD);
310:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[z(tf)]/d[y0]  d[z(tf)]/d[z0]\n");
311:   VecView(user.lambda[1],PETSC_VIEWER_STDOUT_WORLD);
312:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameters: d[y(tf)]/d[mu]\n");
313:   VecView(user.mup[0],PETSC_VIEWER_STDOUT_WORLD);
314:   PetscPrintf(PETSC_COMM_WORLD,"\n sensivitity wrt parameters: d[z(tf)]/d[mu]\n");
315:   VecView(user.mup[1],PETSC_VIEWER_STDOUT_WORLD);

317:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318:      Free work space.  All PETSc objects should be destroyed when they
319:      are no longer needed.
320:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
321:   MatDestroy(&user.A);
322:   MatDestroy(&user.Jacp);
323:   VecDestroy(&user.x);
324:   VecDestroy(&user.lambda[0]);
325:   VecDestroy(&user.lambda[1]);
326:   VecDestroy(&user.mup[0]);
327:   VecDestroy(&user.mup[1]);
328:   TSDestroy(&ts);

330:   PetscFinalize();
331:   return(0);
332: }