Actual source code: ex16opt_p.c

petsc-3.6.1 2015-08-06
Report Typos and Errors
  2: static char help[] = "Solves the van der Pol equation.\n\
  3: Input parameters include:\n\
  4:       -mu : stiffness parameter\n\n";

  6: /*
  7:    Concepts: TS^time-dependent nonlinear problems
  8:    Concepts: TS^van der Pol equation
  9:    Concepts: Optimization using adjoint sensitivity analysis
 10:    Processors: 1
 11: */
 12: /* ------------------------------------------------------------------------

 14:    This program solves the van der Pol equation
 15:        y'' - \mu (1-y^2)*y' + y = 0        (1)
 16:    on the domain 0 <= x <= 1, with the boundary conditions
 17:        y(0) = 2, y'(0) = 0,
 18:    This is a nonlinear equation.

 20:    Notes:
 21:    This code demonstrates the TS solver interface to two variants of
 22:    linear problems, u_t = f(u,t), namely turning (1) into a system of
 23:    first order differential equations,

 25:    [ y' ] = [          z          ]
 26:    [ z' ]   [ \mu (1 - y^2) z - y ]

 28:    which then we can write as a vector equation

 30:    [ u_1' ] = [             u_2           ]  (2)
 31:    [ u_2' ]   [ \mu (1 - u_1^2) u_2 - u_1 ]

 33:    which is now in the desired form of u_t = f(u,t). 
 34:   ------------------------------------------------------------------------- */
 35: #include <petsctao.h>
 36: #include <petscts.h>
 37: #include <petscmat.h>
 38: typedef struct _n_User *User;
 39: struct _n_User {
 40:   PetscReal mu;
 41:   PetscReal next_output;
 42:   PetscInt  steps;
 43:   PetscReal ftime,x_ob[2];
 44:   Mat       A;             /* Jacobian matrix */
 45:   Mat       Jacp;          /* JacobianP matrix */
 46:   Vec       x,lambda[2],mup[2];        /* adjoint variables */
 47: };

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

 51: /*
 52: *  User-defined routines
 53: */
 56: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec X,Vec F,void *ctx)
 57: {
 58:   PetscErrorCode    ierr;
 59:   User              user = (User)ctx;
 60:   PetscScalar       *f;
 61:   const PetscScalar *x;

 64:   VecGetArrayRead(X,&x);
 65:   VecGetArray(F,&f);
 66:   f[0] = x[1];
 67:   f[1] = user->mu*(1.-x[0]*x[0])*x[1]-x[0];
 68:   VecRestoreArrayRead(X,&x);
 69:   VecRestoreArray(F,&f);
 70:   return(0);
 71: }

 75: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx)
 76: {
 77:   PetscErrorCode    ierr;
 78:   User              user = (User)ctx;
 79:   PetscReal         mu   = user->mu;
 80:   PetscInt          rowcol[] = {0,1};
 81:   PetscScalar       J[2][2];
 82:   const PetscScalar *x;

 85:   VecGetArrayRead(X,&x);
 86:   J[0][0] = 0;
 87:   J[1][0] = -2.*mu*x[1]*x[0]-1.;
 88:   J[0][1] = 1.0;
 89:   J[1][1] = mu*(1.0-x[0]*x[0]);
 90:   MatSetValues(A,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 91:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 92:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 93:   if (A != B) {
 94:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 95:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
 96:   }
 97:   VecRestoreArrayRead(X,&x);
 98:   return(0);
 99: }

103: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx)
104: {
105:   PetscErrorCode    ierr;
106:   PetscInt          row[] = {0,1},col[]={0};
107:   PetscScalar       J[2][1];
108:   const PetscScalar *x;

111:   VecGetArrayRead(X,&x);
112:   J[0][0] = 0;
113:   J[1][0] = (1.-x[0]*x[0])*x[1];
114:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
115:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
116:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
117:   VecRestoreArrayRead(X,&x);
118:   return(0);
119: }

123: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
124: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
125: {
126:   PetscErrorCode    ierr;
127:   const PetscScalar *x;
128:   PetscReal         tfinal, dt, tprev;
129:   User              user = (User)ctx;

132:   TSGetTimeStep(ts,&dt);
133:   TSGetDuration(ts,NULL,&tfinal);
134:   TSGetPrevTime(ts,&tprev);
135:   VecGetArrayRead(X,&x);
136:   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]));
137:   PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
138:   VecRestoreArrayRead(X,&x);
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:   PetscMPIInt        size;
151:   struct _n_User     user;
152:   PetscErrorCode     ierr;
153:   Tao                tao;
154:   TaoConvergedReason reason;
155:   Vec                lowerb,upperb;
156:   KSP                ksp;
157:   PC                 pc;

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

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

167:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
168:     Set runtime options
169:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
170:   user.mu          = 1.0;
171:   user.next_output = 0.0;
172:   user.steps       = 0;
173:   user.ftime       = 0.5;

175:   PetscOptionsGetReal(NULL,"-mu",&user.mu,NULL);
176:   PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);

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

187:   MatCreate(PETSC_COMM_WORLD,&user.Jacp);
188:   MatSetSizes(user.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
189:   MatSetFromOptions(user.Jacp);
190:   MatSetUp(user.Jacp);
191: 
192:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193:      Create timestepping solver context
194:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195:   TSCreate(PETSC_COMM_WORLD,&ts);
196:   TSSetType(ts,TSRK);
197:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
198:   TSSetDuration(ts,PETSC_DEFAULT,user.ftime);
199:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
200:   if (monitor) {
201:     TSMonitorSet(ts,Monitor,&user,NULL);
202:   }

204:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
205:      Set initial conditions
206:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
207:   VecGetArray(user.x,&x_ptr);
208:   x_ptr[0] = 2.0;   x_ptr[1] = 0.66666654321;
209:   VecRestoreArray(user.x,&x_ptr);
210:   TSSetTime(ts,0.0);
211:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime));

213:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214:     Save trajectory of solution so that TSAdjointSolve() may be used
215:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216:   TSSetSaveTrajectory(ts);
217: 
218:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
219:      Set runtime options
220:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
221:   TSSetFromOptions(ts);

223:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
224:      Solve nonlinear system
225:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
226:   TSSolve(ts,user.x);
227:   TSGetSolveTime(ts,&(user.ftime));
228:   TSGetTimeStepNumber(ts,&user.steps);
229:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);

231:   VecGetArray(user.x,&x_ptr);
232:   user.x_ob[0] = x_ptr[0];
233:   user.x_ob[1] = x_ptr[1];

235:   MatCreateVecs(user.A,&user.lambda[0],NULL);
236:   MatCreateVecs(user.A,&user.lambda[1],NULL);
237:   MatCreateVecs(user.Jacp,&user.mup[0],NULL);
238:   MatCreateVecs(user.Jacp,&user.mup[1],NULL);

240:   /* Create TAO solver and set desired solution method */
241:   TaoCreate(PETSC_COMM_WORLD,&tao);
242:   TaoSetType(tao,TAOCG);

244:   /* Set initial solution guess */
245:   MatCreateVecs(user.Jacp,&p,NULL);
246:   VecGetArray(p,&x_ptr);
247:   x_ptr[0]   = 6.0;
248:   VecRestoreArray(p,&x_ptr);

250:   TaoSetInitialVector(tao,p);

252:   /* Set routine for function and gradient evaluation */
253:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);
254: 
255:   VecDuplicate(p,&lowerb);
256:   VecDuplicate(p,&upperb);
257:   VecGetArray(lowerb,&x_ptr);
258:   x_ptr[0] = 0.0;
259:   VecRestoreArray(lowerb,&x_ptr);
260:   VecGetArray(upperb,&x_ptr);
261:   x_ptr[0] = 100.0;
262:   VecRestoreArray(upperb,&x_ptr);

264:   TaoSetVariableBounds(tao,lowerb,upperb);
265: 
266:   /* Check for any TAO command line options */
267:   TaoSetFromOptions(tao);
268:   TaoGetKSP(tao,&ksp);
269:   if (ksp) {
270:     KSPGetPC(ksp,&pc);
271:     PCSetType(pc,PCNONE);
272:   }
273: 
274:   TaoSetTolerances(tao,1e-13,1e-13,1e-13,PETSC_DEFAULT,PETSC_DEFAULT);

276:   /* SOLVE THE APPLICATION */
277:   TaoSolve(tao);

279:   /* Get information on termination */
280:   TaoGetConvergedReason(tao,&reason);
281:   if (reason <= 0){
282:       ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
283:   }

285:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
286:   /* Free TAO data structures */
287:   TaoDestroy(&tao);

289:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
290:      Free work space.  All PETSc objects should be destroyed when they
291:      are no longer needed.
292:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
293:   MatDestroy(&user.A);
294:   MatDestroy(&user.Jacp);
295:   VecDestroy(&user.x);
296:   VecDestroy(&user.lambda[0]);
297:   VecDestroy(&user.lambda[1]);
298:   VecDestroy(&user.mup[0]);
299:   VecDestroy(&user.mup[1]);
300:   TSDestroy(&ts);

302:   VecDestroy(&lowerb);
303:   VecDestroy(&upperb);
304:   VecDestroy(&p);
305:   PetscFinalize();
306:   return(0);
307: }

309: /* ------------------------------------------------------------------ */
312: /*
313:    FormFunctionGradient - Evaluates the function and corresponding gradient.

315:    Input Parameters:
316:    tao - the Tao context
317:    X   - the input vector
318:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

320:    Output Parameters:
321:    f   - the newly evaluated function
322:    G   - the newly evaluated gradient
323: */
324: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx)
325: {
326:   User              user = (User)ctx;
327:   TS                ts;
328:   PetscScalar       *x_ptr,*y_ptr;
329:   PetscErrorCode    ierr;

331:   VecGetArray(P,&x_ptr);
332:   user->mu = x_ptr[0];

334:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
335:      Create timestepping solver context
336:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
337:   TSCreate(PETSC_COMM_WORLD,&ts);
338:   TSSetType(ts,TSRK);
339:   TSSetRHSFunction(ts,NULL,RHSFunction,user);
340:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

342:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
343:      Set initial conditions
344:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
345:   VecGetArray(user->x,&x_ptr);
346:   x_ptr[0] = 2;   x_ptr[1] = 0.66666654321;
347:   VecRestoreArray(user->x,&x_ptr);
348:   TSSetTime(ts,0.0);
349:   TSSetDuration(ts,PETSC_DEFAULT,0.5);

351:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
352:     Save trajectory of solution so that TSAdjointSolve() may be used
353:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
354:   TSSetSaveTrajectory(ts);

356:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
357:      Set runtime options
358:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
359:   TSSetFromOptions(ts);

361:   TSSolve(ts,user->x);
362:   TSGetSolveTime(ts,&user->ftime);
363:   TSGetTimeStepNumber(ts,&user->steps);
364:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime);

366:   VecGetArray(user->x,&x_ptr);
367:   *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]);
368:   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));

370:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
371:      Adjoint model starts here
372:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
373:   /*   Redet initial conditions for the adjoint integration */
374:   VecGetArray(user->lambda[0],&y_ptr);
375:   y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
376:   y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
377:   VecRestoreArray(user->lambda[0],&y_ptr);
378:   VecGetArray(user->lambda[1],&x_ptr);
379:   x_ptr[0] = 0.0;   x_ptr[1] = 1.0;
380:   VecRestoreArray(user->lambda[1],&x_ptr);

382:   VecGetArray(user->mup[0],&x_ptr);
383:   x_ptr[0] = 0.0;
384:   VecRestoreArray(user->mup[0],&x_ptr);
385:   VecGetArray(user->mup[1],&x_ptr);
386:   x_ptr[0] = 0.0;
387:   VecRestoreArray(user->mup[1],&x_ptr);
388:   TSSetCostGradients(ts,1,user->lambda,user->mup);

390:   TSSetRHSJacobian(ts,user->A,user->A,RHSJacobian,user);
391:   TSAdjointSetRHSJacobian(ts,user->Jacp,RHSJacobianP,user);

393:   TSAdjointSolve(ts);

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

397:   TSDestroy(&ts);
398:   return(0);
399: }