Actual source code: ex16opt_ic.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: static char help[] = "Solves an ODE-constrained optimization problem -- finding the optimal initial conditions 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 optimal values for initial conditions.
 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;

 26:   PetscInt  steps;
 27:   PetscReal ftime,x_ob[2];
 28:   Mat       A;             /* Jacobian matrix */
 29:   Vec       x,lambda[2];   /* adjoint variables */
 30: };

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

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

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

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

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

 86: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
 87: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
 88: {
 89:   PetscErrorCode    ierr;
 90:   const PetscScalar *x;
 91:   PetscReal         tfinal, dt, tprev;
 92:   User              user = (User)ctx;

 95:   TSGetTimeStep(ts,&dt);
 96:   TSGetDuration(ts,NULL,&tfinal);
 97:   TSGetPrevTime(ts,&tprev);
 98:   VecGetArrayRead(X,&x);
 99:   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]));
100:   PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
101:   VecGetArrayRead(X,&x);
102:   return(0);
103: }

107: int main(int argc,char **argv)
108: {
109:   TS                 ts;          /* nonlinear solver */
110:   Vec                ic;
111:   PetscBool          monitor = PETSC_FALSE;
112:   PetscScalar        *x_ptr;
113:   PetscMPIInt        size;
114:   struct _n_User     user;
115:   PetscErrorCode     ierr;
116:   Tao                tao;
117:   KSP                ksp;
118:   PC                 pc;

120:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
121:      Initialize program
122:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
123:   PetscInitialize(&argc,&argv,NULL,help);

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

128:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
129:     Set runtime options
130:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
131:   user.mu          = 1.0;
132:   user.next_output = 0.0;
133:   user.steps       = 0;
134:   user.ftime       = 0.5;

136:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);
137:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);

139:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
140:     Create necessary matrix and vectors, solve same ODE on every process
141:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
142:   MatCreate(PETSC_COMM_WORLD,&user.A);
143:   MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
144:   MatSetFromOptions(user.A);
145:   MatSetUp(user.A);
146:   MatCreateVecs(user.A,&user.x,NULL);

148:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149:      Create timestepping solver context
150:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151:   TSCreate(PETSC_COMM_WORLD,&ts);
152:   TSSetType(ts,TSRK);
153:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
154:   TSSetDuration(ts,PETSC_DEFAULT,user.ftime);
155:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
156:   if (monitor) {
157:     TSMonitorSet(ts,Monitor,&user,NULL);
158:   }

160:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
161:      Set initial conditions
162:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
163:   VecGetArray(user.x,&x_ptr);
164:   x_ptr[0] = 2.0;   x_ptr[1] = 0.66666654321;
165:   VecRestoreArray(user.x,&x_ptr);
166:   TSSetTime(ts,0.0);
167:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime));

169:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
170:     Save trajectory of solution so that TSAdjointSolve() may be used
171:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
172:   TSSetSaveTrajectory(ts);

174:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
175:      Set runtime options
176:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
177:   TSSetFromOptions(ts);

179:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180:      Solve nonlinear system
181:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182:   TSSolve(ts,user.x);
183:   TSGetSolveTime(ts,&(user.ftime));
184:   TSGetTimeStepNumber(ts,&user.steps);
185:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);

187:   VecGetArray(user.x,&x_ptr);
188:   user.x_ob[0] = x_ptr[0];
189:   user.x_ob[1] = x_ptr[1];

191:   MatCreateVecs(user.A,&user.lambda[0],NULL);

193:   /* Create TAO solver and set desired solution method */
194:   TaoCreate(PETSC_COMM_WORLD,&tao);
195:   TaoSetType(tao,TAOCG);

197:   /* Set initial solution guess */
198:   MatCreateVecs(user.A,&ic,NULL);
199:   VecGetArray(ic,&x_ptr);
200:   x_ptr[0]  = 2.1;
201:   x_ptr[1]  = 0.7;
202:   VecRestoreArray(ic,&x_ptr);

204:   TaoSetInitialVector(tao,ic);

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

209:   /* Check for any TAO command line options */
210:   TaoSetFromOptions(tao);
211:   TaoGetKSP(tao,&ksp);
212:   if (ksp) {
213:     KSPGetPC(ksp,&pc);
214:     PCSetType(pc,PCNONE);
215:   }

217:   TaoSetTolerances(tao,1e-10,PETSC_DEFAULT,PETSC_DEFAULT);

219:   /* SOLVE THE APPLICATION */
220:   TaoSolve(tao);

222:   /* Free TAO data structures */
223:   TaoDestroy(&tao);

225:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226:      Free work space.  All PETSc objects should be destroyed when they
227:      are no longer needed.
228:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229:   MatDestroy(&user.A);
230:   VecDestroy(&user.x);
231:   VecDestroy(&user.lambda[0]);
232:   TSDestroy(&ts);

234:   VecDestroy(&ic);
235:   PetscFinalize();
236:   return(0);
237: }

239: /* ------------------------------------------------------------------ */
242: /*
243:    FormFunctionGradient - Evaluates the function and corresponding gradient.

245:    Input Parameters:
246:    tao - the Tao context
247:    X   - the input vector
248:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

250:    Output Parameters:
251:    f   - the newly evaluated function
252:    G   - the newly evaluated gradient
253: */
254: PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
255: {
256:   User              user = (User)ctx;
257:   TS                ts;
258:   PetscScalar       *x_ptr,*y_ptr;
259:   PetscErrorCode    ierr;
260:   PetscScalar       *ic_ptr;

262:   VecCopy(IC,user->x);

264:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265:      Create timestepping solver context
266:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
267:   TSCreate(PETSC_COMM_WORLD,&ts);
268:   TSSetType(ts,TSRK);
269:   TSSetRHSFunction(ts,NULL,RHSFunction,user);

271:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
272:      Set time
273:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
274:   TSSetTime(ts,0.0);
275:   TSSetInitialTimeStep(ts,0.0,.001);
276:   TSSetDuration(ts,2000,0.5);
277:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

279:   TSSetTolerances(ts,1e-7,NULL,1e-7,NULL);

281:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
282:     Save trajectory of solution so that TSAdjointSolve() may be used
283:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
284:   TSSetSaveTrajectory(ts);

286:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
287:      Set runtime options
288:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
289:   TSSetFromOptions(ts);

291:   TSSolve(ts,user->x);
292:   TSGetSolveTime(ts,&user->ftime);
293:   TSGetTimeStepNumber(ts,&user->steps);
294:   PetscPrintf(PETSC_COMM_WORLD,"mu %.6f, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime);

296:   VecGetArray(IC,&ic_ptr);
297:   VecGetArray(user->x,&x_ptr);
298:   *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]);
299:   PetscPrintf(PETSC_COMM_WORLD,"Observed value y_ob=[%f; %f], ODE solution y=[%f;%f], Cost function f=%f\n",(double)user->x_ob[0],(double)user->x_ob[1],(double)x_ptr[0],(double)x_ptr[1],(double)(*f));

301:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
302:      Adjoint model starts here
303:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
304:   /*   Redet initial conditions for the adjoint integration */
305:   VecGetArray(user->lambda[0],&y_ptr);
306:   y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
307:   y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
308:   VecRestoreArray(user->lambda[0],&y_ptr);
309:   TSSetCostGradients(ts,1,user->lambda,NULL);

311:   /*   Set RHS Jacobian  for the adjoint integration */
312:   TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user);

314:   TSAdjointSolve(ts);

316:   VecCopy(user->lambda[0],G);

318:   TSDestroy(&ts);
319:   return(0);
320: }