Actual source code: ex16opt_ic.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;

 43:   PetscInt  steps;
 44:   PetscReal ftime,x_ob[2];
 45:   Mat       A;             /* Jacobian matrix */
 46:   Vec       x,lambda[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: /* Monitor timesteps and use interpolation to output at integer multiples of 0.1 */
104: static PetscErrorCode Monitor(TS ts,PetscInt step,PetscReal t,Vec X,void *ctx)
105: {
106:   PetscErrorCode    ierr;
107:   const PetscScalar *x;
108:   PetscReal         tfinal, dt, tprev;
109:   User              user = (User)ctx;

112:   TSGetTimeStep(ts,&dt);
113:   TSGetDuration(ts,NULL,&tfinal);
114:   TSGetPrevTime(ts,&tprev);
115:   VecGetArrayRead(X,&x);
116:   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]));
117:   PetscPrintf(PETSC_COMM_WORLD,"t %.6f (tprev = %.6f) \n",(double)t,(double)tprev);
118:   VecGetArrayRead(X,&x);
119:   return(0);
120: }

124: int main(int argc,char **argv)
125: {
126:   TS                 ts;          /* nonlinear solver */
127:   Vec                ic;
128:   PetscBool          monitor = PETSC_FALSE;
129:   PetscScalar        *x_ptr;
130:   PetscMPIInt        size;
131:   struct _n_User     user;
132:   PetscErrorCode     ierr;
133:   Tao                tao;
134:   TaoConvergedReason reason;
135:   KSP                ksp;
136:   PC                 pc;

138:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
139:      Initialize program
140:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
141:   PetscInitialize(&argc,&argv,NULL,help);

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

146:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
147:     Set runtime options
148:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
149:   user.mu          = 1.0;
150:   user.next_output = 0.0;
151:   user.steps       = 0;
152:   user.ftime       = 0.5;

154:   PetscOptionsGetReal(NULL,"-mu",&user.mu,NULL);
155:   PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);

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

166:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167:      Create timestepping solver context
168:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169:   TSCreate(PETSC_COMM_WORLD,&ts);
170:   TSSetType(ts,TSRK);
171:   TSSetRHSFunction(ts,NULL,RHSFunction,&user);
172:   TSSetDuration(ts,PETSC_DEFAULT,user.ftime);
173:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
174:   if (monitor) {
175:     TSMonitorSet(ts,Monitor,&user,NULL);
176:   }

178:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
179:      Set initial conditions
180:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
181:   VecGetArray(user.x,&x_ptr);
182:   x_ptr[0] = 2.0;   x_ptr[1] = 0.66666654321;
183:   VecRestoreArray(user.x,&x_ptr);
184:   TSSetTime(ts,0.0);
185:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)(user.ftime));

187:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
188:     Save trajectory of solution so that TSAdjointSolve() may be used
189:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
190:   TSSetSaveTrajectory(ts);

192:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
193:      Set runtime options
194:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
195:   TSSetFromOptions(ts);

197:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
198:      Solve nonlinear system
199:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
200:   TSSolve(ts,user.x);
201:   TSGetSolveTime(ts,&(user.ftime));
202:   TSGetTimeStepNumber(ts,&user.steps);
203:   PetscPrintf(PETSC_COMM_WORLD,"mu %g, steps %D, ftime %g\n",(double)user.mu,user.steps,(double)user.ftime);

205:   VecGetArray(user.x,&x_ptr);
206:   user.x_ob[0] = x_ptr[0];
207:   user.x_ob[1] = x_ptr[1];

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

211:   /* Create TAO solver and set desired solution method */
212:   TaoCreate(PETSC_COMM_WORLD,&tao);
213:   TaoSetType(tao,TAOCG);

215:   /* Set initial solution guess */
216:   MatCreateVecs(user.A,&ic,NULL);
217:   VecGetArray(ic,&x_ptr);
218:   x_ptr[0]  = 2.1;
219:   x_ptr[1]  = 0.7;
220:   VecRestoreArray(ic,&x_ptr);
221: 
222:   TaoSetInitialVector(tao,ic);

224:   /* Set routine for function and gradient evaluation */
225:   TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&user);
226: 
227:   /* Check for any TAO command line options */
228:   TaoSetFromOptions(tao);
229:   TaoGetKSP(tao,&ksp);
230:   if (ksp) {
231:     KSPGetPC(ksp,&pc);
232:     PCSetType(pc,PCNONE);
233:   }
234: 
235:   TaoSetTolerances(tao,1e-10,1e-10,1e-10,PETSC_DEFAULT,PETSC_DEFAULT);

237:   /* SOLVE THE APPLICATION */
238:   TaoSolve(tao);

240:   /* Get information on termination */
241:   TaoGetConvergedReason(tao,&reason);
242:   if (reason <= 0){
243:       ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
244:   }
245: 
246:   /* Free TAO data structures */
247:   TaoDestroy(&tao);

249:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
250:      Free work space.  All PETSc objects should be destroyed when they
251:      are no longer needed.
252:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
253:   MatDestroy(&user.A);
254:   VecDestroy(&user.x);
255:   VecDestroy(&user.lambda[0]);
256:   TSDestroy(&ts);

258:   VecDestroy(&ic);
259:   PetscFinalize();
260:   return(0);
261: }

263: /* ------------------------------------------------------------------ */
266: /*
267:    FormFunctionGradient - Evaluates the function and corresponding gradient.

269:    Input Parameters:
270:    tao - the Tao context
271:    X   - the input vector
272:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

274:    Output Parameters:
275:    f   - the newly evaluated function
276:    G   - the newly evaluated gradient
277: */
278: PetscErrorCode FormFunctionGradient(Tao tao,Vec IC,PetscReal *f,Vec G,void *ctx)
279: {
280:   User              user = (User)ctx;
281:   TS                ts;
282:   PetscScalar       *x_ptr,*y_ptr;
283:   PetscErrorCode    ierr;
284:   PetscScalar       *ic_ptr;

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

288:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
289:      Create timestepping solver context
290:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
291:   TSCreate(PETSC_COMM_WORLD,&ts);
292:   TSSetType(ts,TSRK);
293:   TSSetRHSFunction(ts,NULL,RHSFunction,user);
294:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
295: 
296:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
297:      Set time
298:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
299:   TSSetTime(ts,0.0);
300:   TSSetInitialTimeStep(ts,0.0,.001);
301:   TSSetDuration(ts,2000,0.5);

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

305:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
306:     Save trajectory of solution so that TSAdjointSolve() may be used
307:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
308:   TSSetSaveTrajectory(ts);
309: 
310:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
311:      Set runtime options
312:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
313:   TSSetFromOptions(ts);

315:   TSSolve(ts,user->x);
316:   TSGetSolveTime(ts,&user->ftime);
317:   TSGetTimeStepNumber(ts,&user->steps);
318:   PetscPrintf(PETSC_COMM_WORLD,"mu %.6f, steps %D, ftime %g\n",(double)user->mu,user->steps,(double)user->ftime);

320:   VecGetArray(IC,&ic_ptr);
321:   VecGetArray(user->x,&x_ptr);
322:   *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]);
323:   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));

325:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
326:      Adjoint model starts here
327:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
328:   /*   Redet initial conditions for the adjoint integration */
329:   VecGetArray(user->lambda[0],&y_ptr);
330:   y_ptr[0] = 2.*(x_ptr[0]-user->x_ob[0]);
331:   y_ptr[1] = 2.*(x_ptr[1]-user->x_ob[1]);
332:   VecRestoreArray(user->lambda[0],&y_ptr);
333:   TSSetCostGradients(ts,1,user->lambda,NULL);

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

338:   TSAdjointSolve(ts);

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

342:   TSDestroy(&ts);
343:   return(0);
344: }