Actual source code: ex16adj.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: static char help[] = "Performs adjoint sensitivity analysis for the van der Pol equation.\n\
  2: Input parameters include:\n\
  3:       -mu : stiffness parameter\n\n";

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

 13:    This program solves the van der Pol equation
 14:        y'' - \mu (1-y^2)*y' + y = 0        (1)
 15:    on the domain 0 <= x <= 1, with the boundary conditions
 16:        y(0) = 2, y'(0) = 0,
 17:    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 adjoint.

 19:    Notes:
 20:    This code demonstrates the TSAdjoint interface to a system of ordinary differential equations (ODEs) in the form of u_t = F(u,t).

 22:    (1) can be turned into a system of first order ODEs
 23:    [ y' ] = [          z          ]
 24:    [ z' ]   [ \mu (1 - y^2) z - y ]

 26:    which then we can write as a vector equation

 28:    [ u_1' ] = [             u_2           ]  (2)
 29:    [ u_2' ]   [ \mu (1 - u_1^2) u_2 - u_1 ]

 31:    which is now in the form of u_t = F(u,t).

 33:    The user provides the right-hand-side function

 35:    [ F(u,t) ] = [ u_2                       ]
 36:                 [ \mu (1 - u_1^2) u_2 - u_1 ]

 38:    the Jacobian function

 40:    dF   [       0           ;         1        ]
 41:    -- = [                                      ]
 42:    du   [ -2 \mu u_1*u_2 - 1;  \mu (1 - u_1^2) ]

 44:    and the JacobainP (the Jacobian w.r.t. parameter) function

 46:    dF      [       0          ]
 47:    ---   = [                  ]
 48:    d\mu    [ (1 - u_1^2) u_2  ]


 51:   ------------------------------------------------------------------------- */

 53: #include <petscts.h>
 54: #include <petscmat.h>
 55: typedef struct _n_User *User;
 56: struct _n_User {
 57:   PetscReal mu;
 58:   PetscReal next_output;
 59:   PetscReal tprev;
 60: };

 62: /*
 63: *  User-defined routines
 64: */
 67: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
 68: {
 69:   PetscErrorCode    ierr;
 70:   User              user = (User)ctx;
 71:   PetscScalar       *f;
 72:   const PetscScalar *x;

 75:   VecGetArrayRead(X,&x);
 76:   VecGetArray(F,&f);
 77:   f[0] = x[1];
 78:   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
 79:   VecRestoreArrayRead(X,&x);
 80:   VecRestoreArray(F,&f);
 81:   return(0);
 82: }

 86: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
 87: {
 88:   PetscErrorCode    ierr;
 89:   User              user = (User)ctx;
 90:   PetscReal         mu   = user->mu;
 91:   PetscInt          rowcol[] = {0,1};
 92:   PetscScalar       J[2][2];
 93:   const PetscScalar *x;

 96:   VecGetArrayRead(X,&x);
 97:   J[0][0] = 0;
 98:   J[1][0] = -2.*mu*x[1]*x[0]-1.;
 99:   J[0][1] = 1.0;
100:   J[1][1] = mu*(1.0-x[0]*x[0]);
101:   MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
102:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
103:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
104:   if (A != B) {
105:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
106:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
107:   }
108:   VecRestoreArrayRead(X,&x);
109:   return(0);
110: }

114: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
115: {
116:   PetscErrorCode    ierr;
117:   PetscInt          row[] = {0,1},col[]={0};
118:   PetscScalar       J[2][1];
119:   const PetscScalar *x;

122:   VecGetArrayRead(X,&x);
123:   J[0][0] = 0;
124:   J[1][0] = (1.-x[0]*x[0])*x[1];
125:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
126:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
127:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
128:   VecRestoreArrayRead(X,&x);
129:   return(0);
130: }

134: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
135: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
136: {
137:   PetscErrorCode    ierr;
138:   const PetscScalar *x;
139:   PetscReal         tfinal, dt, tprev;
140:   User              user = (User)ctx;

143:   TSGetTimeStep(ts,&dt);
144:   TSGetDuration(ts,NULL,&tfinal);
145:   TSGetPrevTime(ts,&tprev);
146:   VecGetArrayRead(X,&x);
147:   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]));
148:   PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
149:   VecRestoreArrayRead(X,&x);
150:   return(0);
151: }

155: int main(int argc,char **argv)
156: {
157:   TS             ts;            /* nonlinear solver */
158:   Vec            x;             /* solution, residual vectors */
159:   Mat            A;             /* Jacobian matrix */
160:   Mat            Jacp;          /* JacobianP matrix */
161:   PetscInt       steps;
162:   PetscReal      ftime   =0.5;
163:   PetscBool      monitor = PETSC_FALSE;
164:   PetscScalar    *x_ptr;
165:   PetscMPIInt    size;
166:   struct _n_User user;
168:   Vec            lambda[2],mu[2];

170:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171:      Initialize program
172:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173:   PetscInitialize(&argc,&argv,NULL,help);

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

178:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179:     Set runtime options
180:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181:   user.mu          = 1;
182:   user.next_output = 0.0;


185:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
186:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);

188:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
189:     Create necessary matrix and vectors, solve same ODE on every process
190:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
191:   MatCreate(PETSC_COMM_WORLD,&A);
192:   MatSetSizes(A,PETSC_DECIDE,PETSC_DECIDE,2,2);
193:   MatSetFromOptions(A);
194:   MatSetUp(A);
195:   MatCreateVecs(A,&x,NULL);

197:   MatCreate(PETSC_COMM_WORLD,&Jacp);
198:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
199:   MatSetFromOptions(Jacp);
200:   MatSetUp(Jacp);

202:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203:      Create timestepping solver context
204:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205:   TSCreate(PETSC_COMM_WORLD,&ts);
206:   TSSetType(ts,TSRK);
207:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
208:   TSSetDuration(ts,PETSC_DEFAULT,ftime);
209:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
210:   if (monitor) {
211:     TSMonitorSet(ts,Monitor,&user,NULL);
212:   }

214:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215:      Set initial conditions
216:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217:   VecGetArray(x,&x_ptr);

219:   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
220:   VecRestoreArray(x,&x_ptr);
221:   TSSetInitialTimeStep(ts,0.0,.001);

223:   /*
224:     Have the TS save its trajectory so that TSAdjointSolve() may be used
225:   */
226:   TSSetSaveTrajectory(ts);

228:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229:      Set runtime options
230:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231:   TSSetFromOptions(ts);

233:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
234:      Solve nonlinear system
235:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236:   TSSolve(ts,x);
237:   TSGetSolveTime(ts,&ftime);
238:   TSGetTimeStepNumber(ts,&steps);
239:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,steps,(double)ftime);
240:   VecView(x,PETSC_VIEWER_STDOUT_WORLD);

242:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
243:      Start the Adjoint model
244:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
245:   MatCreateVecs(A,&lambda[0],NULL);
246:   MatCreateVecs(A,&lambda[1],NULL);
247:   /*   Reset initial conditions for the adjoint integration */
248:   VecGetArray(lambda[0],&x_ptr);
249:   x_ptr[0] = 1.0;   x_ptr[1] = 0.0;
250:   VecRestoreArray(lambda[0],&x_ptr);
251:   VecGetArray(lambda[1],&x_ptr);
252:   x_ptr[0] = 0.0;   x_ptr[1] = 1.0;
253:   VecRestoreArray(lambda[1],&x_ptr);

255:   MatCreateVecs(Jacp,&mu[0],NULL);
256:   MatCreateVecs(Jacp,&mu[1],NULL);
257:   VecGetArray(mu[0],&x_ptr);
258:   x_ptr[0] = 0.0;
259:   VecRestoreArray(mu[0],&x_ptr);
260:   VecGetArray(mu[1],&x_ptr);
261:   x_ptr[0] = 0.0;
262:   VecRestoreArray(mu[1],&x_ptr);
263:   TSSetCostGradients(ts,2,lambda,mu);

265:   /*   Set RHS Jacobian for the adjoint integration */
266:   TSSetRHSJacobian(ts,A,A,RHSJacobian,&user);

268:   /*   Set RHS JacobianP */
269:   TSAdjointSetRHSJacobian(ts,Jacp,RHSJacobianP,&user);

271:   TSAdjointSolve(ts);

273:   VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
274:   VecView(lambda[1],PETSC_VIEWER_STDOUT_WORLD);
275:   VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
276:   VecView(mu[1],PETSC_VIEWER_STDOUT_WORLD);

278:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
279:      Free work space.  All PETSc objects should be destroyed when they
280:      are no longer needed.
281:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
282:   MatDestroy(&A);
283:   MatDestroy(&Jacp);
284:   VecDestroy(&x);
285:   VecDestroy(&lambda[0]);
286:   VecDestroy(&lambda[1]);
287:   VecDestroy(&mu[0]);
288:   VecDestroy(&mu[1]);
289:   TSDestroy(&ts);

291:   PetscFinalize();
292:   return(0);
293: }