Actual source code: ex20opt_ic.c

petsc-3.7.7 2017-09-25
Report Typos and Errors
  1: #define c11 1.0
  2: #define c12 0
  3: #define c21 2.0
  4: #define c22 1.0
  5: static char help[] = "Solves a DAE-constrained optimization problem -- finding the optimal initial conditions for the van der Pol equation.\n";

  7: /*
  8:    Concepts: TS^time-dependent nonlinear problems
  9:    Concepts: TS^van der Pol equation DAE equivalent
 10:    Concepts: Optimization using adjoint sensitivities
 11:    Processors: 1
 12: */
 13: /* ------------------------------------------------------------------------
 14:    Notes:
 15:    This code demonstrates how to solve a DAE-constrained optimization problem with TAO, TSAdjoint and TS.
 16:    The nonlinear problem is written in a DAE equivalent form.
 17:    The objective is to minimize the difference between observation and model prediction by finding optimal values for initial conditions.
 18:    The gradient is computed with the discrete adjoint of an implicit theta method, see ex20adj.c for details.
 19:   ------------------------------------------------------------------------- */
 20: #include <petsctao.h>
 21: #include <petscts.h>

 23: typedef struct _n_User *User;
 24: struct _n_User {
 25:   PetscReal mu;
 26:   PetscReal next_output;

 28:   /* Sensitivity analysis support */
 29:   PetscReal ftime,x_ob[2];
 30:   Mat       A;            /* Jacobian matrix */
 31:   Vec       x,lambda[2];  /* adjoint variables */
 32: };

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

 36: /*
 37: *  User-defined routines
 38: */
 41: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 42: {
 43:   PetscErrorCode    ierr;
 44:   User              user = (User)ctx;
 45:   PetscScalar       *f;
 46:   const PetscScalar *x,*xdot;

 49:   VecGetArrayRead(X,&x);
 50:   VecGetArrayRead(Xdot,&xdot);
 51:   VecGetArray(F,&f);
 52:   f[0] = xdot[0] - x[1];
 53:   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
 54:   VecRestoreArrayRead(X,&x);
 55:   VecRestoreArrayRead(Xdot,&xdot);
 56:   VecRestoreArray(F,&f);
 57:   return(0);
 58: }

 62: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 63: {
 64:   PetscErrorCode    ierr;
 65:   User              user     = (User)ctx;
 66:   PetscInt          rowcol[] = {0,1};
 67:   PetscScalar       J[2][2];
 68:   const PetscScalar *x;

 71:   VecGetArrayRead(X,&x);

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

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

 79:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 80:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 81:   if (A != B) {
 82:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 83:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 84:   }
 85:   return(0);
 86: }

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

100:   TSGetTimeStep(ts,&dt);
101:   TSGetDuration(ts,NULL,&tfinal);

103:   while (user->next_output <= t && user->next_output <= tfinal) {
104:     VecDuplicate(X,&interpolatedX);
105:     TSInterpolate(ts,user->next_output,interpolatedX);
106:     VecGetArrayRead(interpolatedX,&x);
107:     PetscPrintf(PETSC_COMM_WORLD,"[%.1f] %D TS %.6f (dt = %.6f) X % 12.6e % 12.6e\n",
108:                        user->next_output,step,t,dt,(double)PetscRealPart(x[0]),
109:                        (double)PetscRealPart(x[1]));
110:     VecRestoreArrayRead(interpolatedX,&x);
111:     VecDestroy(&interpolatedX);
112:     user->next_output += 0.1;
113:   }
114:   return(0);
115: }

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

132:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
133:      Initialize program
134:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
135:   PetscInitialize(&argc,&argv,NULL,help);

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

140:   /* Create TAO solver and set desired solution method */
141:   TaoCreate(PETSC_COMM_WORLD,&tao);
142:   TaoSetType(tao,TAOCG);

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:     Set runtime options
146:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147:   user.next_output = 0.0;
148:   user.mu          = 1.0e6;
149:   user.ftime       = 0.5;
150:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
151:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);

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

162:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
163:      Create timestepping solver context
164:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
165:   TSCreate(PETSC_COMM_WORLD,&ts);
166:   TSSetType(ts,TSCN);
167:   TSSetIFunction(ts,NULL,IFunction,&user);
168:   TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);
169:   TSSetDuration(ts,PETSC_DEFAULT,user.ftime);
170:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

172:  if (monitor) {
173:     TSMonitorSet(ts,Monitor,&user,NULL);
174:   }

176:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
177:      Set initial conditions
178:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
179:   VecGetArray(user.x,&x_ptr);
180:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
181:   VecRestoreArray(user.x,&x_ptr);
182:   TSSetInitialTimeStep(ts,0.0,.0001);

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);

199:   VecGetArray(user.x,&x_ptr);
200:   user.x_ob[0] = x_ptr[0];
201:   user.x_ob[1] = x_ptr[1];

203:   /* Create sensitivity variable */
204:   MatCreateVecs(user.A,&user.lambda[0],NULL);
205:   MatCreateVecs(user.A,&user.lambda[1],NULL);

207:   /* Set initial solution guess */
208:   MatCreateVecs(user.A,&ic,NULL);
209:   VecGetArray(ic,&x_ptr);
210:   x_ptr[0] = 2.1;
211:   x_ptr[1] = -0.66666654321;
212:   VecRestoreArray(ic,&x_ptr);

214:   TaoSetInitialVector(tao,ic);

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

219:   /* Check for any TAO command line options */
220:   TaoSetFromOptions(tao);
221:   TaoGetKSP(tao,&ksp);
222:   if (ksp) {
223:     KSPGetPC(ksp,&pc);
224:     PCSetType(pc,PCNONE);
225:   }

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

229:   /* SOLVE THE APPLICATION */
230:   TaoSolve(tao);

232:   VecView(ic,PETSC_VIEWER_STDOUT_WORLD);
233:   /* Free TAO data structures */
234:   TaoDestroy(&tao);

236:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237:      Free work space.  All PETSc objects should be destroyed when they
238:      are no longer needed.
239:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240:   MatDestroy(&user.A);
241:   VecDestroy(&user.x);
242:   VecDestroy(&user.lambda[0]);
243:   VecDestroy(&user.lambda[1]);
244:   TSDestroy(&ts);

246:   VecDestroy(&ic);
247:   PetscFinalize();
248:   return(0);
249: }


252: /* ------------------------------------------------------------------ */
255: /*
256:    FormFunctionGradient - Evaluates the function and corresponding gradient.

258:    Input Parameters:
259:    tao - the Tao context
260:    X   - the input vector
261:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

263:    Output Parameters:
264:    f   - the newly evaluated function
265:    G   - the newly evaluated gradient
266: */
267: PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
268: {
269:   User           user_ptr = (User)ctx;
270:   TS             ts;
271:   PetscScalar    *x_ptr,*y_ptr;

274:   VecCopy(IC,user_ptr->x);

276:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
277:      Create timestepping solver context
278:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
279:   TSCreate(PETSC_COMM_WORLD,&ts);
280:   TSSetType(ts,TSCN);
281:   TSSetIFunction(ts,NULL,IFunction,user_ptr);
282:   TSSetIJacobian(ts,user_ptr->A,user_ptr->A,IJacobian,user_ptr);

284:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
285:      Set time
286:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
287:   TSSetTime(ts,0.0);
288:   TSSetDuration(ts,PETSC_DEFAULT,0.5);
289:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

291:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
292:     Save trajectory of solution so that TSAdjointSolve() may be used
293:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
294:   TSSetSaveTrajectory(ts);

296:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
297:      Set runtime options
298:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
299:   TSSetFromOptions(ts);

301:   TSSolve(ts,user_ptr->x);
302:   VecGetArray(user_ptr->x,&x_ptr);
303:   *f   = (x_ptr[0]-user_ptr->x_ob[0])*(x_ptr[0]-user_ptr->x_ob[0])+(x_ptr[1]-user_ptr->x_ob[1])*(x_ptr[1]-user_ptr->x_ob[1]);
304:   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)x_ptr[0],(double)x_ptr[1],(double)(*f));

306:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
307:      Adjoint model starts here
308:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
309:   /*   Redet initial conditions for the adjoint integration */
310:   VecGetArray(user_ptr->lambda[0],&y_ptr);
311:   y_ptr[0] = 2.*(x_ptr[0]-user_ptr->x_ob[0]);
312:   y_ptr[1] = 2.*(x_ptr[1]-user_ptr->x_ob[1]);
313:   VecRestoreArray(user_ptr->lambda[0],&y_ptr);
314:   TSSetCostGradients(ts,1,user_ptr->lambda,NULL);

316:   TSAdjointSolve(ts);
317:   VecCopy(user_ptr->lambda[0],G);
318:   TSDestroy(&ts);
319:   return(0);
320: }