Actual source code: ex16opt_p.c

petsc-3.11.4 2019-09-28
Report Typos and Errors
  1: static char help[] = "Solves an ODE-constrained optimization problem -- finding the optimal stiffness parameter 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: Optimization using adjoint sensitivities
  9:    Processors: 1
 10: */
 11: /* ------------------------------------------------------------------------

 13:    Notes:
 14:    This code demonstrates how to solve an ODE-constrained optimization problem with TAO, TSAdjoint and TS.
 15:    The objective is to minimize the difference between observation and model prediction by finding an optimal value for parameter \mu.
 16:    The gradient is computed with the discrete adjoint of an explicit Runge-Kutta method, see ex16adj.c for details.
 17:   ------------------------------------------------------------------------- */
 18:  #include <petsctao.h>
 19:  #include <petscts.h>
 20:  #include <petscmat.h>
 21: typedef struct _n_User *User;
 22: struct _n_User {
 23:   PetscReal mu;
 24:   PetscReal next_output;
 25:   PetscInt  steps;
 26:   PetscReal ftime,x_ob[2];
 27:   Mat       A;                  /* Jacobian matrix */
 28:   Mat       Jacp;               /* JacobianP matrix */
 29:   Vec       x,lambda[2],mup[2]; /* adjoint variables */
 30: };

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

 34: /*
 35: *  User-defined routines
 36: */
 37: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
 38: {
 39:   PetscErrorCode    ierr;
 40:   User              user = (User)ctx;
 41:   PetscScalar       *f;
 42:   const PetscScalar *x;

 45:   VecGetArrayRead(X,&x);
 46:   VecGetArray(F,&f);
 47:   f[0] = x[1];
 48:   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
 49:   VecRestoreArrayRead(X,&x);
 50:   VecRestoreArray(F,&f);
 51:   return(0);
 52: }

 54: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
 55: {
 56:   PetscErrorCode    ierr;
 57:   User              user = (User)ctx;
 58:   PetscReal         mu   = user->mu;
 59:   PetscInt          rowcol[] = {0,1};
 60:   PetscScalar       J[2][2];
 61:   const PetscScalar *x;

 64:   VecGetArrayRead(X,&x);
 65:   J[0][0] = 0;
 66:   J[1][0] = -2.*mu*x[1]*x[0]-1.;
 67:   J[0][1] = 1.0;
 68:   J[1][1] = mu*(1.0-x[0]*x[0]);
 69:   MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 70:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 71:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 72:   if (B && A != B) {
 73:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 74:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 75:   }
 76:   VecRestoreArrayRead(X,&x);
 77:   return(0);
 78: }

 80: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
 81: {
 82:   PetscErrorCode    ierr;
 83:   PetscInt          row[] = {0,1},col[]={0};
 84:   PetscScalar       J[2][1];
 85:   const PetscScalar *x;

 88:   VecGetArrayRead(X,&x);
 89:   J[0][0] = 0;
 90:   J[1][0] = (1.-x[0]*x[0])*x[1];
 91:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
 92:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 93:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 94:   VecRestoreArrayRead(X,&x);
 95:   return(0);
 96: }

 98: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
 99: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
100: {
101:   PetscErrorCode    ierr;
102:   const PetscScalar *x;
103:   PetscReal         tfinal, dt, tprev;
104:   User              user = (User)ctx;

107:   TSGetTimeStep(ts,&dt);
108:   TSGetMaxTime(ts,&tfinal);
109:   TSGetPrevTime(ts,&tprev);
110:   VecGetArrayRead(X,&x);
111:   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]));
112:   PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
113:   VecRestoreArrayRead(X,&x);
114:   return(0);
115: }

117: int main(int argc,char **argv)
118: {
119:   TS                 ts;          /* nonlinear solver */
120:   Vec                p;
121:   PetscBool          monitor = PETSC_FALSE;
122:   PetscScalar        *x_ptr;
123:   PetscMPIInt        size;
124:   struct _n_User     user;
125:   PetscErrorCode     ierr;
126:   Tao                tao;
127:   Vec                lowerb,upperb;
128:   KSP                ksp;
129:   PC                 pc;

131:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
132:      Initialize program
133:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
134:   PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
135:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
136:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

138:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139:     Set runtime options
140:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141:   user.mu          = 1.0;
142:   user.next_output = 0.0;
143:   user.steps       = 0;
144:   user.ftime       = 0.5;

146:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
147:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);

149:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
150:     Create necessary matrix and vectors, solve same ODE on every process
151:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
152:   MatCreate(PETSC_COMM_WORLD,&user.A);
153:   MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
154:   MatSetFromOptions(user.A);
155:   MatSetUp(user.A);
156:   MatCreateVecs(user.A,&user.x,NULL);

158:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
159:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
160:   MatSetFromOptions(user.Jacp);
161:   MatSetUp(user.Jacp);

163:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
164:      Create timestepping solver context
165:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
166:   TSCreate(PETSC_COMM_WORLD,&ts);
167:   TSSetType(ts,TSRK);
168:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
169:   TSSetMaxTime(ts,user.ftime);
170:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
171:   if (monitor) {
172:     TSMonitorSet(ts,Monitor,&user,NULL);
173:   }

175:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176:      Set initial conditions
177:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178:   VecGetArray(user.x,&x_ptr);
179:   x_ptr[0] = 2.0;   x_ptr[1] = 0.66666654321;
180:   VecRestoreArray(user.x,&x_ptr);
181:   TSSetTime(ts,0.0);
182:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime));

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185:     Save trajectory of solution so that TSAdjointSolve() may be used
186:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187:   TSSetSaveTrajectory(ts);

189:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
190:      Set runtime options
191:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
192:   TSSetFromOptions(ts);

194:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195:      Solve nonlinear system
196:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197:   TSSolve(ts,user.x);
198:   TSGetSolveTime(ts,&(user.ftime));
199:   TSGetStepNumber(ts,&user.steps);
200:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);
201:   TSDestroy(&ts);

203:   VecGetArray(user.x,&x_ptr);
204:   user.x_ob[0] = x_ptr[0];
205:   user.x_ob[1] = x_ptr[1];
206:   VecRestoreArray(user.x,&x_ptr);

208:   MatCreateVecs(user.A,&user.lambda[0],NULL);
209:   MatCreateVecs(user.A,&user.lambda[1],NULL);
210:   MatCreateVecs(user.Jacp,&user.mup[0],NULL);
211:   MatCreateVecs(user.Jacp,&user.mup[1],NULL);

213:   /* Create TAO solver and set desired solution method */
214:   TaoCreate(PETSC_COMM_WORLD,&tao);
215:   TaoSetType(tao,TAOCG);

217:   /* Set initial solution guess */
218:   MatCreateVecs(user.Jacp,&p,NULL);
219:   VecGetArray(p,&x_ptr);
220:   x_ptr[0]   = 6.0;
221:   VecRestoreArray(p,&x_ptr);

223:   TaoSetInitialVector(tao,p);

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

228:   VecDuplicate(p,&lowerb);
229:   VecDuplicate(p,&upperb);
230:   VecGetArray(lowerb,&x_ptr);
231:   x_ptr[0] = 0.0;
232:   VecRestoreArray(lowerb,&x_ptr);
233:   VecGetArray(upperb,&x_ptr);
234:   x_ptr[0] = 100.0;
235:   VecRestoreArray(upperb,&x_ptr);

237:   TaoSetVariableBounds(tao,lowerb,upperb);

239:   /* Check for any TAO command line options */
240:   TaoSetFromOptions(tao);
241:   TaoGetKSP(tao,&ksp);
242:   if (ksp) {
243:     KSPGetPC(ksp,&pc);
244:     PCSetType(pc,PCNONE);
245:   }

247:   TaoSetTolerances(tao,1e-13,PETSC_DEFAULT,PETSC_DEFAULT);

249:   /* SOLVE THE APPLICATION */
250:   TaoSolve(tao);

252:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
253:   /* Free TAO data structures */
254:   TaoDestroy(&tao);

256:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
257:      Free work space.  All PETSc objects should be destroyed when they
258:      are no longer needed.
259:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
260:   MatDestroy(&user.A);
261:   MatDestroy(&user.Jacp);
262:   VecDestroy(&user.x);
263:   VecDestroy(&user.lambda[0]);
264:   VecDestroy(&user.lambda[1]);
265:   VecDestroy(&user.mup[0]);
266:   VecDestroy(&user.mup[1]);

268:   VecDestroy(&lowerb);
269:   VecDestroy(&upperb);
270:   VecDestroy(&p);
271:   PetscFinalize();
272:   return ierr;
273: }

275: /* ------------------------------------------------------------------ */
276: /*
277:    FormFunctionGradient - Evaluates the function and corresponding gradient.

279:    Input Parameters:
280:    tao - the Tao context
281:    X   - the input vector
282:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

284:    Output Parameters:
285:    f   - the newly evaluated function
286:    G   - the newly evaluated gradient
287: */
288: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
289: {
290:   User              user = (User)ctx;
291:   TS                ts;
292:   PetscScalar       *x_ptr,*y_ptr;
293:   PetscErrorCode    ierr;

296:   VecGetArrayRead(P,(const PetscScalar**)&x_ptr);
297:   user->mu = x_ptr[0];
298:   VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr);

300:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
301:      Create timestepping solver context
302:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
303:   TSCreate(PETSC_COMM_WORLD,&ts);
304:   TSSetType(ts,TSRK);
305:   TSSetRHSFunction(ts,NULL,RHSFunction,user);

307:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
308:      Set initial conditions
309:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
310:   VecGetArray(user->x,&x_ptr);
311:   x_ptr[0] = 2;
312:   x_ptr[1] = 0.66666654321;
313:   VecRestoreArray(user->x,&x_ptr);
314:   TSSetTime(ts,0.0);
315:   TSSetMaxTime(ts,0.5);
316:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

318:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
319:     Save trajectory of solution so that TSAdjointSolve() may be used
320:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
321:   TSSetSaveTrajectory(ts);

323:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
324:      Set runtime options
325:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
326:   TSSetFromOptions(ts);

328:   TSSolve(ts,user->x);
329:   TSGetSolveTime(ts,&user->ftime);
330:   TSGetStepNumber(ts,&user->steps);
331:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime);

333:   VecGetArray(user->x,&x_ptr);
334:   *f   = (x_ptr[0]-user->x_ob[0])*(x_ptr[0]-user->x_ob[0])+(x_ptr[1]-user->x_ob[1])*(x_ptr[1]-user->x_ob[1]);
335:   PetscPrintf(PETSC_COMM_WORLD,"Observed value y_ob=[%f; %f], ODE solution y=%f, Cost function f=%f\n",(double)user->x_ob[0],(double)user->x_ob[1],(double)x_ptr[0],(double)(*f));

337:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
338:      Adjoint model starts here
339:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
340:   /*   Redet initial conditions for the adjoint integration */
341:   VecGetArray(user->lambda[0],&y_ptr);
342:   y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
343:   y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
344:   VecRestoreArray(user->x,&x_ptr);
345:   VecRestoreArray(user->lambda[0],&y_ptr);
346:   VecRestoreArray(user->x,&x_ptr);
347:   VecGetArray(user->lambda[1],&x_ptr);
348:   x_ptr[0] = 0.0;
349:   x_ptr[1] = 1.0;
350:   VecRestoreArray(user->lambda[1],&x_ptr);
351:   VecGetArray(user->mup[0],&x_ptr);
352:   x_ptr[0] = 0.0;
353:   VecRestoreArray(user->mup[0],&x_ptr);
354:   VecGetArray(user->mup[1],&x_ptr);
355:   x_ptr[0] = 0.0;
356:   VecRestoreArray(user->mup[1],&x_ptr);
357:   TSSetCostGradients(ts,1,user->lambda,user->mup);

359:   TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user);
360:   TSSetRHSJacobianP(ts,user->Jacp,RHSJacobianP,user);

362:   TSAdjointSolve(ts);

364:   VecCopy(user->mup[0],G);

366:   TSDestroy(&ts);
367:   return(0);
368: }

370: /*TEST

372:     build:
373:       requires: !complex !single

375:     test:
376:       suffix: 1
377:       args: -monitor 0 -viewer_binary_skip_info -tao_view -tao_monitor  -tao_gttol 1.e-5 -ts_trajectory_dirname ex16opt_pdir

379: TEST*/