Actual source code: ex20opt_ic.c

petsc-3.6.4 2016-04-12
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 the van der Pol equation.\n\
  6: Input parameters include:\n";

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

 16:    This program solves the van der Pol DAE ODE equivalent
 17:        y' = z                 (1)
 18:        z' = mu[(1-y^2)z-y]
 19:    on the domain 0 <= x <= 1, with the boundary conditions
 20:        y(0) = 2, y'(0) = -6.666665432100101e-01,
 21:    and
 22:        mu = 10^6.
 23:    This is a nonlinear equation.

 25:    Notes:
 26:    This code demonstrates the TS solver interface to a variant of
 27:    linear problems, u_t = f(u,t), namely turning (1) into a system of
 28:    first order differential equations,

 30:    [ y' ] = [          z          ]
 31:    [ z' ]   [     mu[(1-y^2)z-y]  ]

 33:    which then we can write as a vector equation

 35:    [ u_1' ] = [      u_2              ]  (2)
 36:    [ u_2' ]   [ mu[(1-u_1^2)u_2-u_1]  ]

 38:    which is now in the desired form of u_t = f(u,t). 

 40:   ------------------------------------------------------------------------- */
 41: #include <petsctao.h>
 42: #include <petscts.h>

 44: typedef struct _n_User *User;
 45: struct _n_User {
 46:   PetscReal mu;
 47:   PetscReal next_output;
 48: 
 49:   /* Sensitivity analysis support */
 50:   PetscReal ftime,x_ob[2];
 51:   Mat       A;            /* Jacobian matrix */
 52:   Vec       x,lambda[2];  /* adjoint variables */
 53: };

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

 57: /*
 58: *  User-defined routines
 59: */
 62: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,void *ctx)
 63: {
 64:   PetscErrorCode    ierr;
 65:   User              user = (User)ctx;
 66:   PetscScalar       *f;
 67:   const PetscScalar *x,*xdot;

 70:   VecGetArrayRead(X,&x);
 71:   VecGetArrayRead(Xdot,&xdot);
 72:   VecGetArray(F,&f);
 73:   f[0] = xdot[0] - x[1];
 74:   f[1] = c21*(xdot[0]-x[1]) + xdot[1] - user->mu*((1.0-x[0]*x[0])*x[1] - x[0]) ;
 75:   VecRestoreArrayRead(X,&x);
 76:   VecRestoreArrayRead(Xdot,&xdot);
 77:   VecRestoreArray(F,&f);
 78:   return(0);
 79: }

 83: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,void *ctx)
 84: {
 85:   PetscErrorCode    ierr;
 86:   User              user     = (User)ctx;
 87:   PetscInt          rowcol[] = {0,1};
 88:   PetscScalar       J[2][2];
 89:   const PetscScalar *x;

 92:   VecGetArrayRead(X,&x);

 94:   J[0][0] = a;     J[0][1] =  -1.0;
 95:   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]);
 96: 
 97:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 98:   VecRestoreArrayRead(X,&x);

100:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
101:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
102:   if (A != B) {
103:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
104:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
105:   }
106:   return(0);
107: }

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

121:   TSGetTimeStep(ts,&dt);
122:   TSGetDuration(ts,NULL,&tfinal);

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

140: int main(int argc,char **argv)
141: {
142:   TS                 ts;            /* nonlinear solver */
143:   Vec                ic;
144:   PetscBool          monitor = PETSC_FALSE;
145:   PetscScalar        *x_ptr;
146:   PetscMPIInt        size;
147:   struct _n_User     user;
148:   PetscErrorCode     ierr;
149:   Tao                tao;
150:   TaoConvergedReason reason;
151:   KSP                ksp;
152:   PC                 pc;

154:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
155:      Initialize program
156:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
157:   PetscInitialize(&argc,&argv,NULL,help);
158: 
159:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
160:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

162:   /* Create TAO solver and set desired solution method */
163:   TaoCreate(PETSC_COMM_WORLD,&tao);
164:   TaoSetType(tao,TAOCG);

166:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
167:     Set runtime options
168:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
169:   user.next_output = 0.0;
170:   user.mu          = 1.0e6;
171:   user.ftime       = 0.5;
172:   PetscOptionsGetBool(NULL,"-monitor",&monitor,NULL);
173:   PetscOptionsGetReal(NULL,"-mu",&user.mu,NULL);

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

184:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
185:      Create timestepping solver context
186:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
187:   TSCreate(PETSC_COMM_WORLD,&ts);
188:   TSSetType(ts,TSCN);
189:   TSSetIFunction(ts,NULL,IFunction,&user);
190:   TSSetIJacobian(ts,user.A,user.A,IJacobian,&user);
191:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
192:   TSSetDuration(ts,PETSC_DEFAULT,user.ftime);
193: 
194:  if (monitor) {
195:     TSMonitorSet(ts,Monitor,&user,NULL);
196:   }

198:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
199:      Set initial conditions
200:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
201:   VecGetArray(user.x,&x_ptr);
202:   x_ptr[0] = 2.0;   x_ptr[1] = -0.66666654321;
203:   VecRestoreArray(user.x,&x_ptr);
204:   TSSetInitialTimeStep(ts,0.0,.0001);

206:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
207:     Save trajectory of solution so that TSAdjointSolve() may be used
208:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
209:   TSSetSaveTrajectory(ts);

211:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
212:      Set runtime options
213:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
214:   TSSetFromOptions(ts);

216:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
217:      Solve nonlinear system
218:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
219:   TSSolve(ts,user.x);

221:   VecGetArray(user.x,&x_ptr);
222:   user.x_ob[0] = x_ptr[0];
223:   user.x_ob[1] = x_ptr[1];

225:   /* Create sensitivity variable */
226:   MatCreateVecs(user.A,&user.lambda[0],NULL);
227:   MatCreateVecs(user.A,&user.lambda[1],NULL);

229:   /* Set initial solution guess */
230:   MatCreateVecs(user.A,&ic,NULL);
231:   VecGetArray(ic,&x_ptr);
232:   x_ptr[0] = 2.1;
233:   x_ptr[1] = -0.66666654321;
234:   VecRestoreArray(ic,&x_ptr);

236:   TaoSetInitialVector(tao,ic);

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

241:   /* Check for any TAO command line options */
242:   TaoSetFromOptions(tao);
243:   TaoGetKSP(tao,&ksp);
244:   if (ksp) {
245:     KSPGetPC(ksp,&pc);
246:     PCSetType(pc,PCNONE);
247:   }

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

251:   /* SOLVE THE APPLICATION */
252:   TaoSolve(tao);

254:   /* Get information on termination */
255:   TaoGetConvergedReason(tao,&reason);
256:   if (reason <= 0){
257:       ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
258:   }

260:   VecView(ic,PETSC_VIEWER_STDOUT_WORLD);
261:   /* Free TAO data structures */
262:   TaoDestroy(&tao);

264:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
265:      Free work space.  All PETSc objects should be destroyed when they
266:      are no longer needed.
267:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
268:   MatDestroy(&user.A);
269:   VecDestroy(&user.x);
270:   VecDestroy(&user.lambda[0]);
271:   VecDestroy(&user.lambda[1]);
272:   TSDestroy(&ts);

274:   VecDestroy(&ic);
275:   PetscFinalize();
276:   return(0);
277: }


280: /* ------------------------------------------------------------------ */
283: /*
284:    FormFunctionGradient - Evaluates the function and corresponding gradient.

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

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

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

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

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

319:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
320:     Save trajectory of solution so that TSAdjointSolve() may be used
321:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
322:   TSSetSaveTrajectory(ts);
323: 
324:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
325:      Set runtime options
326:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
327:   TSSetFromOptions(ts);

329:   TSSolve(ts,user_ptr->x);
330:   VecGetArray(user_ptr->x,&x_ptr);
331:   *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]);
332:   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));

334:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
335:      Adjoint model starts here
336:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
337:   /*   Redet initial conditions for the adjoint integration */
338:   VecGetArray(user_ptr->lambda[0],&y_ptr);
339:   y_ptr[0] = 2.*(x_ptr[0]-user_ptr->x_ob[0]);
340:   y_ptr[1] = 2.*(x_ptr[1]-user_ptr->x_ob[1]);
341:   VecRestoreArray(user_ptr->lambda[0],&y_ptr);
342:   TSSetCostGradients(ts,1,user_ptr->lambda,NULL);

344:   TSAdjointSolve(ts);
345:   VecCopy(user_ptr->lambda[0],G);
346:   TSDestroy(&ts);
347:   return(0);
348: }