Actual source code: ex20opt_ic.c

petsc-3.10.5 2019-03-28
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: */
 39: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 40: {
 41:   PetscErrorCode    ierr;
 42:   User              user = (User)ctx;
 43:   PetscScalar       *f;
 44:   const PetscScalar *x,*xdot;

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

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

 67:   VecGetArrayRead(X,&x);

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

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

 75:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 76:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 77:   if (A != B) {
 78:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 79:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 80:   }
 81:   return(0);
 82: }

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

 94:   TSGetTimeStep(ts,&dt);
 95:   TSGetMaxTime(ts,&tfinal);

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

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

124:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
125:      Initialize program
126:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
127:   PetscInitialize(&argc,&argv,NULL,help);
128:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
129:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

131:   /* Create TAO solver and set desired solution method */
132:   TaoCreate(PETSC_COMM_WORLD,&tao);
133:   TaoSetType(tao,TAOBLMVM);

135:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
136:     Set runtime options
137:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
138:   user.next_output = 0.0;
139:   user.mu          = 1.0;
140:   user.ftime       = 1.0;
141:   PetscOptionsGetBool(NULL,NULL,"-monitor",&monitor,NULL);
142:   PetscOptionsGetReal(NULL,NULL,"-mu",&user.mu,NULL);

144:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
145:     Create necessary matrix and vectors, solve same ODE on every process
146:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
147:   MatCreate(PETSC_COMM_WORLD,&user.A);
148:   MatSetSizes(user.A,PETSC_DECIDE,PETSC_DECIDE,2,2);
149:   MatSetFromOptions(user.A);
150:   MatSetUp(user.A);
151:   MatCreateVecs(user.A,&user.x,NULL);

153:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
154:      Create timestepping solver context
155:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
156:   TSCreate(PETSC_COMM_WORLD,&ts);
157:   TSSetType(ts,TSCN);
158:   TSSetIFunction(ts,NULL,IFunction,&user);
159:   TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);
160:   TSSetMaxTime(ts,user.ftime);
161:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

163:  if (monitor) {
164:     TSMonitorSet(ts,Monitor,&user,NULL);
165:   }

167:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168:      Set initial conditions
169:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170:   VecGetArray(user.x,&x_ptr);
171:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
172:   VecRestoreArray(user.x,&x_ptr);
173:   TSSetTimeStep(ts,0.03125);

175:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
176:     Save trajectory of solution so that TSAdjointSolve() may be used
177:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
178:   TSSetSaveTrajectory(ts);

180:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181:      Set runtime options
182:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183:   TSSetFromOptions(ts);
184:   TSSolve(ts,user.x);

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

190:   /* Create sensitivity variable */
191:   MatCreateVecs(user.A,&user.lambda[0],NULL);
192:   MatCreateVecs(user.A,&user.lambda[1],NULL);

194:   /* Set initial solution guess */
195:   MatCreateVecs(user.A,&ic,NULL);
196:   VecGetArray(ic,&x_ptr);
197:   x_ptr[0] = 2.2;
198:   x_ptr[1] = -0.7;
199:   VecRestoreArray(ic,&x_ptr);

201:   TaoSetInitialVector(tao,ic);

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

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

214:   /* SOLVE THE APPLICATION */
215:   TaoSolve(tao);

217:   VecView(ic,PETSC_VIEWER_STDOUT_WORLD);
218:   /* Free TAO data structures */
219:   TaoDestroy(&tao);

221:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222:      Free work space.  All PETSc objects should be destroyed when they
223:      are no longer needed.
224:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
225:   MatDestroy(&user.A);
226:   VecDestroy(&user.x);
227:   VecDestroy(&user.lambda[0]);
228:   VecDestroy(&user.lambda[1]);
229:   TSDestroy(&ts);
230:   VecDestroy(&ic);
231:   PetscFinalize();
232:   return ierr;
233: }


236: /* ------------------------------------------------------------------ */
237: /*
238:    FormFunctionGradient - Evaluates the function and corresponding gradient.

240:    Input Parameters:
241:    tao - the Tao context
242:    X   - the input vector
243:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

245:    Output Parameters:
246:    f   - the newly evaluated function
247:    G   - the newly evaluated gradient
248: */
249: PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
250: {
251:   User           user_ptr = (User)ctx;
252:   TS             ts;
253:   PetscScalar    *x_ptr,*y_ptr;

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

259:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
260:      Create timestepping solver context
261:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
262:   TSCreate(PETSC_COMM_WORLD,&ts);
263:   TSSetType(ts,TSCN);
264:   TSSetIFunction(ts,NULL,IFunction,user_ptr);
265:   TSSetIJacobian(ts,user_ptr->A,user_ptr->A,IJacobian,user_ptr);

267:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
268:      Set time
269:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
270:   TSSetTime(ts,0.0);
271:   TSSetMaxTime(ts,user_ptr->ftime);
272:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
273:   TSSetTimeStep(ts,0.03125);

275:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
276:     Save trajectory of solution so that TSAdjointSolve() may be used
277:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
278:   TSSetSaveTrajectory(ts);

280:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
281:      Set runtime options
282:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
283:   TSSetFromOptions(ts);

285:   TSSolve(ts,user_ptr->x);
286:   VecGetArray(user_ptr->x,&x_ptr);
287:   *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]);
288:   PetscPrintf(PETSC_COMM_WORLD,"Observed value y_ob=[%g; %g], ODE solution y=[%g;%g], Cost function f=%g\n",(double)user_ptr->x_ob[0],(double)user_ptr->x_ob[1],(double)x_ptr[0],(double)x_ptr[1],(double)(*f));

290:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
291:      Adjoint model starts here
292:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
293:   /*   Redet initial conditions for the adjoint integration */
294:   VecGetArray(user_ptr->lambda[0],&y_ptr);
295:   y_ptr[0] = 2.*(x_ptr[0]-user_ptr->x_ob[0]);
296:   y_ptr[1] = 2.*(x_ptr[1]-user_ptr->x_ob[1]);
297:   VecRestoreArray(user_ptr->lambda[0],&y_ptr);
298:   TSSetCostGradients(ts,1,user_ptr->lambda,NULL);

300:   TSAdjointSolve(ts);
301:   VecCopy(user_ptr->lambda[0],G);
302:   TSDestroy(&ts);
303:   return(0);
304: }

306: /*TEST
307:     build:
308:       requires: !complex !single
309:     test:
310:       args:  -monitor 0 -ts_type theta -ts_theta_endpoint -ts_theta_theta 0.5 -viewer_binary_skip_info -tao_view -mu 1.0 -ts_trajectory_dirname ex20opt_icdir
311:       output_file: output/ex20opt_ic_1.out

313: TEST*/