Actual source code: ex20opt_p.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: #define rescale 10

  7: static char help[] = "Solves the van der Pol equation.\n\
  8: Input parameters include:\n";

 10: /*
 11:    Concepts: TS^time-dependent nonlinear problems
 12:    Concepts: TS^van der Pol equation DAE equivalent
 13:    Concepts: Optimization using adjoint sensitivity analysis
 14:    Processors: 1
 15: */
 16: /* ------------------------------------------------------------------------

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

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

 32:    [ y' ] = [          z          ]
 33:    [ z' ]   [     mu[(1-y^2)z-y]  ]

 35:    which then we can write as a vector equation

 37:    [ u_1' ] = [      u_2              ]  (2)
 38:    [ u_2' ]   [ mu[(1-u_1^2)u_2-u_1]  ]

 40:    which is now in the desired form of u_t = f(u,t). 
 41:   ------------------------------------------------------------------------- */
 42: #include <petsctao.h>
 43: #include <petscts.h>

 45: typedef struct _n_User *User;
 46: struct _n_User {
 47:   PetscReal mu;
 48:   PetscReal next_output;
 49: 
 50:   /* Sensitivity analysis support */
 51:   PetscReal ftime,x_ob[2];
 52:   Mat       A;                       /* Jacobian matrix */
 53:   Mat       Jacp;                    /* JacobianP matrix */
 54:   Vec       x,lambda[2],mup[2];  /* adjoint variables */
 55: };

 57: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);

 59: /*
 60: *  User-defined routines
 61: */
 64: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 65: {
 66:   PetscErrorCode    ierr;
 67:   User              user = (User)ctx;
 68:   PetscScalar       *f;
 69:   const PetscScalar *x,*xdot;

 72:   VecGetArrayRead(X,&x);
 73:   VecGetArrayRead(Xdot,&xdot);
 74:   VecGetArray(F,&f);
 75:   f[0] = xdot[0] - x[1];
 76:   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
 77:   VecRestoreArrayRead(X,&x);
 78:   VecRestoreArrayRead(Xdot,&xdot);
 79:   VecRestoreArray(F,&f);
 80:   return(0);
 81: }

 85: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 86: {
 87:   PetscErrorCode    ierr;
 88:   User              user     = (User)ctx;
 89:   PetscInt          rowcol[] = {0,1};
 90:   PetscScalar       J[2][2];
 91:   const PetscScalar *x;
 93:   VecGetArrayRead(X,&x);

 95:   J[0][0] = a;     J[0][1] =  -1.0;
 96:   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]);
 97: 
 98:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 99:   VecRestoreArrayRead(X,&x);

101:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
102:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
103:   if (A != B) {
104:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
105:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
106:   }
107:   return(0);
108: }

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

121:   J[0][0] = 0;
122:   J[1][0] = (1.-x[0]*x[0])*x[1]-x[0];
123:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
124:   VecRestoreArrayRead(X,&x);

126:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
127:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
128:   return(0);
129: }

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

143:   TSGetTimeStep(ts,&dt);
144:   TSGetDuration(ts,NULL,&tfinal);

146:   while (user->next_output <= t && user->next_output <= tfinal) {
147:     VecDuplicate(X,&interpolatedX);
148:     TSInterpolate(ts,user->next_output,interpolatedX);
149:     VecGetArrayRead(interpolatedX,&x);
150:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
151:                        user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
152:                        (double)PetscRealPart(x[1]));
153:     VecRestoreArrayRead(interpolatedX,&x);
154:     VecDestroy(&interpolatedX);
155:     user->next_output += 0.1;
156:   }
157:   return(0);
158: }

162: int main(int argc,char **argv)
163: {
164:   TS                 ts;            /* nonlinear solver */
165:   Vec                p;
166:   PetscBool          monitor = PETSC_FALSE;
167:   PetscScalar        *x_ptr;
168:   const PetscScalar  *y_ptr;
169:   PetscMPIInt        size;
170:   struct _n_User     user;
171:   PetscErrorCode     ierr;
172:   Tao                tao;
173:   TaoConvergedReason reason;
174:   KSP                ksp;
175:   PC                 pc;

177:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
178:      Initialize program
179:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
180:   PetscInitialize(&argc,&argv,NULL,help);
181: 
182:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
183:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

185:   /* Create TAO solver and set desired solution method */
186:   TaoCreate(PETSC_COMM_WORLD,&tao);
187:   TaoSetType(tao,TAOCG);

189:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190:     Set runtime options
191:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192:   user.next_output = 0.0;
193:   user.mu          = 1.0;
194:   user.ftime       = 0.5;
195:   PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);
196:   PetscOptionsGetReal(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:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
220:   TSSetDuration(ts,PETSC_DEFAULT,user.ftime);
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:      Set runtime options
235:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
236:   TSSetFromOptions(ts);

238:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
239:      Solve nonlinear system
240:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
241:   TSSolve(ts,user.x);

243:   VecGetArrayRead(user.x,&y_ptr);
244:   user.x_ob[0] = y_ptr[0];
245:   user.x_ob[1] = y_ptr[1];
246:   VecRestoreArrayRead(user.x,&y_ptr);

248:   /* Create sensitivity variable */
249:   MatCreateVecs(user.A,&user.lambda[0],NULL);
250:   MatCreateVecs(user.Jacp,&user.mup[0],NULL);

252:   /*
253:      Optimization starts
254:   */
255:   /* Set initial solution guess */
256:   MatCreateVecs(user.Jacp,&p,NULL);
257:   VecGetArray(p,&x_ptr);
258:   x_ptr[0] = 1.2;
259:   VecRestoreArray(p,&x_ptr);
260:   TaoSetInitialVector(tao,p);

262:   /* Set routine for function and gradient evaluation */
263:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);

265:   /* Check for any TAO command line options */
266:   TaoSetFromOptions(tao);
267:   TaoGetKSP(tao,&ksp);
268:   if (ksp) {
269:     KSPGetPC(ksp,&pc);
270:     PCSetType(pc,PCNONE);
271:   }

273:   TaoSetTolerances(tao,1e-7,1e-7,1e-7,1e-7,1e-7);

275:   TaoSolve(tao);
276:   /* Get information on termination */
277:   TaoGetConvergedReason(tao,&reason);
278:   if (reason <= 0){
279:       ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
280:   }
281:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
282:   /* Free TAO data structures */
283:   TaoDestroy(&tao);

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


301: /* ------------------------------------------------------------------ */
304: /*
305:    FormFunctionGradient - Evaluates the function and corresponding gradient.

307:    Input Parameters:
308:    tao - the Tao context
309:    X   - the input vector
310:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

312:    Output Parameters:
313:    f   - the newly evaluated function
314:    G   - the newly evaluated gradient
315: */
316: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
317: {
318:   User              user_ptr = (User)ctx;
319:   TS                ts;
320:   PetscScalar       *x_ptr;
321:   const PetscScalar *y_ptr;
322:   PetscErrorCode    ierr;

324:   VecGetArrayRead(P,&y_ptr);
325:   user_ptr->mu = y_ptr[0];
326:   VecRestoreArrayRead(P,&y_ptr);

328:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
329:      Create timestepping solver context
330:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
331:   TSCreate(PETSC_COMM_WORLD,&ts);
332:   TSSetType(ts,TSCN);
333:   TSSetIFunction(ts,NULL,IFunction,user_ptr);
334:   TSSetIJacobian(ts,user_ptr->A,user_ptr->A,IJacobian,user_ptr);
335:   TSAdjointSetRHSJacobian(ts,user_ptr->Jacp,RHSJacobianP,user_ptr);
336:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

338:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
339:      Set time
340:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
341:   TSSetTime(ts,0.0);
342:   TSSetDuration(ts,PETSC_DEFAULT,0.5);

344:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
345:     Save trajectory of solution so that TSAdjointSolve() may be used
346:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
347:   TSSetSaveTrajectory(ts);

349:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
350:      Set initial conditions
351:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
352:   VecGetArray(user_ptr->x,&x_ptr);
353:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
354:   VecRestoreArray(user_ptr->x,&x_ptr);
355: 
356:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
357:      Set runtime options
358:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
359:   TSSetFromOptions(ts);

361:   TSSolve(ts,user_ptr->x);
362:   VecGetArrayRead(user_ptr->x,&y_ptr);
363:   *f   = rescale*(y_ptr[0]-user_ptr->x_ob[0])*(y_ptr[0]-user_ptr->x_ob[0])+rescale*(y_ptr[1]-user_ptr->x_ob[1])*(y_ptr[1]-user_ptr->x_ob[1]);
364:   PetscPrintf(PETSC_COMM_WORLD,"Observed value y_ob=[%f; %f], ODE solution y=[%f;%f], Cost function f=%f\n",(double)user_ptr->x_ob[0],(double)user_ptr->x_ob[1],(double)y_ptr[0],(double)y_ptr[1],(double)(*f));
365:   /*   Redet initial conditions for the adjoint integration */
366:   VecGetArray(user_ptr->lambda[0],&x_ptr);
367:   x_ptr[0] = rescale*2.*(y_ptr[0]-user_ptr->x_ob[0]);
368:   x_ptr[1] = rescale*2.*(y_ptr[1]-user_ptr->x_ob[1]);
369:   VecRestoreArrayRead(user_ptr->x,&y_ptr);
370:   VecRestoreArray(user_ptr->lambda[0],&x_ptr);

372:   VecGetArray(user_ptr->mup[0],&x_ptr);
373:   x_ptr[0] = 0.0;
374:   VecRestoreArray(user_ptr->mup[0],&x_ptr);
375:   TSSetCostGradients(ts,1,user_ptr->lambda,user_ptr->mup);

377:   TSAdjointSolve(ts);
378:   VecCopy(user_ptr->mup[0],G);
379:   TSDestroy(&ts);
380:   return(0);
381: }