Actual source code: ex20opt_p.c

petsc-3.7.5 2017-01-01
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:   Notes:
 19:   This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
 20:   The nonlinear problem is written in a DAE equivalent form.
 21:   The objective is to minimize the difference between observation and model prediction by finding an optimal value for parameter \mu.
 22:   The gradient is computed with the discrete adjoint of an implicit theta method, see ex20adj.c for details.
 23:   ------------------------------------------------------------------------- */
 24: #include <petsctao.h>
 25: #include <petscts.h>

 27: typedef struct _n_User *User;
 28: struct _n_User {
 29:   PetscReal mu;
 30:   PetscReal next_output;

 32:   /* Sensitivity analysis support */
 33:   PetscReal ftime,x_ob[2];
 34:   Mat       A;                       /* Jacobian matrix */
 35:   Mat       Jacp;                    /* JacobianP matrix */
 36:   Vec       x,lambda[2],mup[2];  /* adjoint variables */
 37: };

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

 41: /*
 42: *  User-defined routines
 43: */
 46: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 47: {
 48:   PetscErrorCode    ierr;
 49:   User              user = (User)ctx;
 50:   PetscScalar       *f;
 51:   const PetscScalar *x,*xdot;

 54:   VecGetArrayRead(X,&x);
 55:   VecGetArrayRead(Xdot,&xdot);
 56:   VecGetArray(F,&f);
 57:   f[0] = xdot[0] - x[1];
 58:   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
 59:   VecRestoreArrayRead(X,&x);
 60:   VecRestoreArrayRead(Xdot,&xdot);
 61:   VecRestoreArray(F,&f);
 62:   return(0);
 63: }

 67: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 68: {
 69:   PetscErrorCode    ierr;
 70:   User              user     = (User)ctx;
 71:   PetscInt          rowcol[] = {0,1};
 72:   PetscScalar       J[2][2];
 73:   const PetscScalar *x;
 75:   VecGetArrayRead(X,&x);

 77:   J[0][0] = a;     J[0][1] =  -1.0;
 78:   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]);

 80:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 81:   VecRestoreArrayRead(X,&x);

 83:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 84:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 85:   if (A != B) {
 86:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 87:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 88:   }
 89:   return(0);
 90: }

 94: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
 95: {
 96:   PetscErrorCode    ierr;
 97:   PetscInt          row[] = {0,1},col[]={0};
 98:   PetscScalar       J[2][1];
 99:   const PetscScalar *x;
101:   VecGetArrayRead(X,&x);

103:   J[0][0] = 0;
104:   J[1][0] = (1.-x[0]*x[0])*x[1]-x[0];
105:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
106:   VecRestoreArrayRead(X,&x);

108:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
109:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
110:   return(0);
111: }

115: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
116: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
117: {
118:   PetscErrorCode    ierr;
119:   const PetscScalar *x;
120:   PetscReal         tfinal, dt;
121:   User              user = (User)ctx;
122:   Vec               interpolatedX;

125:   TSGetTimeStep(ts,&dt);
126:   TSGetDuration(ts,NULL,&tfinal);

128:   while (user->next_output <= t && user->next_output <= tfinal) {
129:     VecDuplicate(X,&interpolatedX);
130:     TSInterpolate(ts,user->next_output,interpolatedX);
131:     VecGetArrayRead(interpolatedX,&x);
132:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
133:                        user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
134:                        (double)PetscRealPart(x[1]));
135:     VecRestoreArrayRead(interpolatedX,&x);
136:     VecDestroy(&interpolatedX);
137:     user->next_output += 0.1;
138:   }
139:   return(0);
140: }

144: int main(int argc,char **argv)
145: {
146:   TS                 ts;            /* nonlinear solver */
147:   Vec                p;
148:   PetscBool          monitor = PETSC_FALSE;
149:   PetscScalar        *x_ptr;
150:   const PetscScalar  *y_ptr;
151:   PetscMPIInt        size;
152:   struct _n_User     user;
153:   PetscErrorCode     ierr;
154:   Tao                tao;
155:   KSP                ksp;
156:   PC                 pc;

158:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
159:      Initialize program
160:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
161:   PetscInitialize(&argc,&argv,NULL,help);

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

166:   /* Create TAO solver and set desired solution method */
167:   TaoCreate(PETSC_COMM_WORLD,&tao);
168:   TaoSetType(tao,TAOCG);

170:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
171:     Set runtime options
172:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
173:   user.next_output = 0.0;
174:   user.mu          = 1.0;
175:   user.ftime       = 0.5;
176:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
177:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);

179:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180:     Create necessary matrix and vectors, solve same ODE on every process
181:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182:   MatCreate(PETSC_COMM_WORLD,&user.A);
183:   MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
184:   MatSetFromOptions(user.A);
185:   MatSetUp(user.A);
186:   MatCreateVecs(user.A,&user.x,NULL);

188:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
189:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
190:   MatSetFromOptions(user.Jacp);
191:   MatSetUp(user.Jacp);

193:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
194:      Create timestepping solver context
195:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
196:   TSCreate(PETSC_COMM_WORLD,&ts);
197:   TSSetType(ts,TSCN);
198:   TSSetIFunction(ts,NULL,IFunction,&user);
199:   TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);
200:   TSSetDuration(ts,PETSC_DEFAULT,user.ftime);
201:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
202:   if (monitor) {
203:     TSMonitorSet(ts,Monitor,&user,NULL);
204:   }

206:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207:      Set initial conditions
208:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209:   VecGetArray(user.x,&x_ptr);
210:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
211:   VecRestoreArray(user.x,&x_ptr);
212:   TSSetInitialTimeStep(ts,0.0,.0001);

214:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
215:      Set runtime options
216:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
217:   TSSetFromOptions(ts);

219:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
220:      Solve nonlinear system
221:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
222:   TSSolve(ts,user.x);

224:   VecGetArrayRead(user.x,&y_ptr);
225:   user.x_ob[0] = y_ptr[0];
226:   user.x_ob[1] = y_ptr[1];
227:   VecRestoreArrayRead(user.x,&y_ptr);

229:   /* Create sensitivity variable */
230:   MatCreateVecs(user.A,&user.lambda[0],NULL);
231:   MatCreateVecs(user.Jacp,&user.mup[0],NULL);

233:   /*
234:      Optimization starts
235:   */
236:   /* Set initial solution guess */
237:   MatCreateVecs(user.Jacp,&p,NULL);
238:   VecGetArray(p,&x_ptr);
239:   x_ptr[0] = 1.2;
240:   VecRestoreArray(p,&x_ptr);
241:   TaoSetInitialVector(tao,p);

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

246:   /* Check for any TAO command line options */
247:   TaoSetFromOptions(tao);
248:   TaoGetKSP(tao,&ksp);
249:   if (ksp) {
250:     KSPGetPC(ksp,&pc);
251:     PCSetType(pc,PCNONE);
252:   }

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

256:   TaoSolve(tao);

258:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
259:   /* Free TAO data structures */
260:   TaoDestroy(&tao);

262:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
263:      Free work space.  All PETSc objects should be destroyed when they
264:      are no longer needed.
265:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
266:   MatDestroy(&user.A);
267:   MatDestroy(&user.Jacp);
268:   VecDestroy(&user.x);
269:   VecDestroy(&user.lambda[0]);
270:   VecDestroy(&user.mup[0]);
271:   TSDestroy(&ts);
272:   VecDestroy(&p);
273:   PetscFinalize();
274:   return(0);
275: }


278: /* ------------------------------------------------------------------ */
281: /*
282:    FormFunctionGradient - Evaluates the function and corresponding gradient.

284:    Input Parameters:
285:    tao - the Tao context
286:    X   - the input vector
287:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

289:    Output Parameters:
290:    f   - the newly evaluated function
291:    G   - the newly evaluated gradient
292: */
293: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
294: {
295:   User              user_ptr = (User)ctx;
296:   TS                ts;
297:   PetscScalar       *x_ptr;
298:   const PetscScalar *y_ptr;
299:   PetscErrorCode    ierr;

301:   VecGetArrayRead(P,&y_ptr);
302:   user_ptr->mu = y_ptr[0];
303:   VecRestoreArrayRead(P,&y_ptr);

305:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
306:      Create timestepping solver context
307:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
308:   TSCreate(PETSC_COMM_WORLD,&ts);
309:   TSSetType(ts,TSCN);
310:   TSSetIFunction(ts,NULL,IFunction,user_ptr);
311:   TSSetIJacobian(ts,user_ptr->A,user_ptr->A,IJacobian,user_ptr);
312:   TSAdjointSetRHSJacobian(ts,user_ptr->Jacp,RHSJacobianP,user_ptr);

314:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
315:      Set time
316:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
317:   TSSetTime(ts,0.0);
318:   TSSetDuration(ts,PETSC_DEFAULT,0.5);
319:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

321:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
322:     Save trajectory of solution so that TSAdjointSolve() may be used
323:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
324:   TSSetSaveTrajectory(ts);

326:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
327:      Set initial conditions
328:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
329:   VecGetArray(user_ptr->x,&x_ptr);
330:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
331:   VecRestoreArray(user_ptr->x,&x_ptr);

333:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
334:      Set runtime options
335:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
336:   TSSetFromOptions(ts);

338:   TSSolve(ts,user_ptr->x);
339:   VecGetArrayRead(user_ptr->x,&y_ptr);
340:   *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]);
341:   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));
342:   /*   Redet initial conditions for the adjoint integration */
343:   VecGetArray(user_ptr->lambda[0],&x_ptr);
344:   x_ptr[0] = rescale*2.*(y_ptr[0]-user_ptr->x_ob[0]);
345:   x_ptr[1] = rescale*2.*(y_ptr[1]-user_ptr->x_ob[1]);
346:   VecRestoreArrayRead(user_ptr->x,&y_ptr);
347:   VecRestoreArray(user_ptr->lambda[0],&x_ptr);

349:   VecGetArray(user_ptr->mup[0],&x_ptr);
350:   x_ptr[0] = 0.0;
351:   VecRestoreArray(user_ptr->mup[0],&x_ptr);
352:   TSSetCostGradients(ts,1,user_ptr->lambda,user_ptr->mup);

354:   TSAdjointSolve(ts);
355:   VecCopy(user_ptr->mup[0],G);
356:   TSDestroy(&ts);
357:   return(0);
358: }