Actual source code: ex3opt_fd.c

petsc-3.7.1 2016-05-15
Report Typos and Errors
  2: static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\n";


\begin{eqnarray}
\frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
\frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
\end{eqnarray}

 13: /*
 14:   Solve the same optimization problem as in ex3opt.c.
 15:   Use finite difference to approximate the gradients.
 16: */
 17: #include <petsctao.h>
 18: #include <petscts.h>

 20: typedef struct {
 21:   PetscScalar H,D,omega_b,omega_s,Pmax,Pmax_ini,Pm,E,V,X,u_s,c;
 22:   PetscInt    beta;
 23:   PetscReal   tf,tcl;
 24: } AppCtx;

 26: PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);

 28: /* Event check */
 31: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
 32: {
 33:   AppCtx        *user=(AppCtx*)ctx;

 36:   /* Event for fault-on time */
 37:   fvalue[0] = t - user->tf;
 38:   /* Event for fault-off time */
 39:   fvalue[1] = t - user->tcl;

 41:   return(0);
 42: }

 46: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
 47: {
 48:   AppCtx *user=(AppCtx*)ctx;


 52:   if (event_list[0] == 0) {
 53:     if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
 54:     else user->Pmax = user->Pmax_ini; /* Going backward, reversal of event */
 55:   } else if(event_list[0] == 1) {
 56:     if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault  - this is done by setting Pmax = Pmax_ini */
 57:     else user->Pmax = 0.0; /* Going backward, reversal of event */
 58:   }
 59:   return(0);
 60: }

 64: /*
 65:      Defines the ODE passed to the ODE solver
 66: */
 67: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 68: {
 69:   PetscErrorCode    ierr;
 70:   PetscScalar       *f,Pmax;
 71:   const PetscScalar *u,*udot;

 74:   /*  The next three lines allow us to access the entries of the vectors directly */
 75:   VecGetArrayRead(U,&u);
 76:   VecGetArrayRead(Udot,&udot);
 77:   VecGetArray(F,&f);
 78:   Pmax = ctx->Pmax;

 80:   f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s);
 81:   f[1] = 2.0*ctx->H/ctx->omega_s*udot[1] +  Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - ctx->omega_s)- ctx->Pm;

 83:   VecRestoreArrayRead(U,&u);
 84:   VecRestoreArrayRead(Udot,&udot);
 85:   VecRestoreArray(F,&f);
 86:   return(0);
 87: }

 91: /*
 92:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 93: */
 94: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
 95: {
 96:   PetscErrorCode    ierr;
 97:   PetscInt          rowcol[] = {0,1};
 98:   PetscScalar       J[2][2],Pmax;
 99:   const PetscScalar *u,*udot;

102:   VecGetArrayRead(U,&u);
103:   VecGetArrayRead(Udot,&udot);
104:   Pmax = ctx->Pmax;

106:   J[0][0] = a;                       J[0][1] = -ctx->omega_b;
107:   J[1][1] = 2.0*ctx->H/ctx->omega_s*a + ctx->D;   J[1][0] = Pmax*PetscCosScalar(u[0]);

109:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
110:   VecRestoreArrayRead(U,&u);
111:   VecRestoreArrayRead(Udot,&udot);

113:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
114:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
115:   if (A != B) {
116:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
117:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
118:   }
119:   return(0);
120: }

124: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
125: {
127:   PetscInt       row[] = {0,1},col[]={0};
128:   PetscScalar    J[2][1];

131:   J[0][0] = 0;
132:   J[1][0] = 1.;
133:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
134:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
135:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
136:   return(0);
137: }

141: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
142: {
143:   PetscErrorCode    ierr;
144:   PetscScalar       *r;
145:   const PetscScalar *u;

148:   VecGetArrayRead(U,&u);
149:   VecGetArray(R,&r);
150:   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
151:   VecRestoreArray(R,&r);
152:   VecRestoreArrayRead(U,&u);
153:   return(0);
154: }

158: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
159: {
160:   PetscErrorCode    ierr;
161:   PetscScalar       *ry;
162:   const PetscScalar *u;

165:   VecGetArrayRead(U,&u);
166:   VecGetArray(drdy[0],&ry);
167:   ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
168:   VecRestoreArray(drdy[0],&ry);
169:   VecRestoreArrayRead(U,&u);
170:   return(0);
171: }

175: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
176: {
178:   PetscScalar    *rp;

181:   VecGetArray(drdp[0],&rp);
182:   rp[0] = 0.;
183:   VecRestoreArray(drdp[0],&rp);
184:   return(0);
185: }

189: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
190: {
191:   PetscErrorCode    ierr;
192:   PetscScalar       *y,sensip;
193:   const PetscScalar *x;

196:   VecGetArrayRead(lambda,&x);
197:   VecGetArray(mu,&y);
198:   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
199:   /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %g \n",(double)sensip); */
200:   y[0] = sensip;
201:   VecRestoreArray(mu,&y);
202:   VecRestoreArrayRead(lambda,&x);
203:   return(0);
204: }

208: int main(int argc,char **argv)
209: {
210:   Vec                p;
211:   PetscScalar        *x_ptr;
212:   PetscErrorCode     ierr;
213:   PetscMPIInt        size;
214:   AppCtx             ctx;
215:   Vec                lowerb,upperb;
216:   Tao                tao;
217:   KSP                ksp;
218:   PC                 pc;

220:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
221:      Initialize program
222:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
223:   PetscInitialize(&argc,&argv,NULL,help);
225:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
226:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

228:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
229:     Set runtime options
230:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
231:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
232:   {
233:     ctx.beta    = 2;
234:     ctx.c       = 10000.0;
235:     ctx.u_s     = 1.0;
236:     ctx.omega_s = 1.0;
237:     ctx.omega_b = 120.0*PETSC_PI;
238:     ctx.H       = 5.0;
239:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
240:     ctx.D       = 5.0;
241:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
242:     ctx.E       = 1.1378;
243:     ctx.V       = 1.0;
244:     ctx.X       = 0.545;
245:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
246:     ctx.Pmax_ini = ctx.Pmax;
247:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
248:     ctx.Pm      = 1.061605;
249:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
250:     ctx.tf      = 0.1;
251:     ctx.tcl     = 0.2;
252:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
253:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);

255:   }
256:   PetscOptionsEnd();

258:   /* Create TAO solver and set desired solution method */
259:   TaoCreate(PETSC_COMM_WORLD,&tao);
260:   TaoSetType(tao,TAOBLMVM);

262:   /*
263:      Optimization starts
264:   */
265:   /* Set initial solution guess */
266:   VecCreateSeq(PETSC_COMM_WORLD,1,&p);
267:   VecGetArray(p,&x_ptr);
268:   x_ptr[0]   = ctx.Pm;
269:   VecRestoreArray(p,&x_ptr);

271:   TaoSetInitialVector(tao,p);
272:   /* Set routine for function and gradient evaluation */
273:   TaoSetObjectiveRoutine(tao,FormFunction,(void *)&ctx);
274:   TaoSetGradientRoutine(tao,TaoDefaultComputeGradient,(void *)&ctx);

276:   /* Set bounds for the optimization */
277:   VecDuplicate(p,&lowerb);
278:   VecDuplicate(p,&upperb);
279:   VecGetArray(lowerb,&x_ptr);
280:   x_ptr[0] = 0.;
281:   VecRestoreArray(lowerb,&x_ptr);
282:   VecGetArray(upperb,&x_ptr);
283:   x_ptr[0] = 1.1;
284:   VecRestoreArray(upperb,&x_ptr);
285:   TaoSetVariableBounds(tao,lowerb,upperb);

287:   /* Check for any TAO command line options */
288:   TaoSetFromOptions(tao);
289:   TaoGetKSP(tao,&ksp);
290:   if (ksp) {
291:     KSPGetPC(ksp,&pc);
292:     PCSetType(pc,PCNONE);
293:   }

295:   TaoSetTolerances(tao,1e-15,1e-15,1e-15);
296:   /* SOLVE THE APPLICATION */
297:   TaoSolve(tao);

299:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);

301:   /* Free TAO data structures */
302:   TaoDestroy(&tao);
303:   VecDestroy(&p);
304:   VecDestroy(&lowerb);
305:   VecDestroy(&upperb);
306:   PetscFinalize();
307:   return 0;
308: }

310: /* ------------------------------------------------------------------ */
313: /*
314:    FormFunction - Evaluates the function and corresponding gradient.

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

321:    Output Parameters:
322:    f   - the newly evaluated function
323: */
324: PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
325: {
326:   AppCtx            *ctx = (AppCtx*)ctx0;
327:   TS                ts;

329:   Vec               U;             /* solution will be stored here */
330:   Mat               A;             /* Jacobian matrix */
331:   Mat               Jacp;          /* Jacobian matrix */
332:   PetscErrorCode    ierr;
333:   PetscInt          n = 2;
334:   PetscReal         ftime;
335:   PetscInt          steps;
336:   PetscScalar       *u;
337:   PetscScalar       *mx_ptr,*y_ptr;
338:   const PetscScalar *x_ptr,*qx_ptr;
339:   Vec               lambda[1],q,mu[1];
340:   PetscInt          direction[2];
341:   PetscBool         terminate[2];

343:   VecGetArrayRead(P,&x_ptr);
344:   ctx->Pm = x_ptr[0];
345:   VecRestoreArrayRead(P,&x_ptr);
346:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
347:     Create necessary matrix and vectors
348:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
349:   MatCreate(PETSC_COMM_WORLD,&A);
350:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
351:   MatSetType(A,MATDENSE);
352:   MatSetFromOptions(A);
353:   MatSetUp(A);

355:   MatCreateVecs(A,&U,NULL);

357:   MatCreate(PETSC_COMM_WORLD,&Jacp);
358:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
359:   MatSetFromOptions(Jacp);
360:   MatSetUp(Jacp);

362:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
363:      Create timestepping solver context
364:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
365:   TSCreate(PETSC_COMM_WORLD,&ts);
366:   TSSetProblemType(ts,TS_NONLINEAR);
367:   TSSetType(ts,TSCN);
368:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
369:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,ctx);
370:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

372:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
373:      Set initial conditions
374:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
375:   VecGetArray(U,&u);
376:   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
377:   u[1] = 1.0;
378:   VecRestoreArray(U,&u);
379:   TSSetSolution(ts,U);

381:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
382:     Save trajectory of solution so that TSAdjointSolve() may be used
383:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
384:   TSSetSaveTrajectory(ts);

386:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
387:      Set solver options
388:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
389:   TSSetDuration(ts,PETSC_DEFAULT,1.0);
390:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_STEPOVER);
391:   TSSetInitialTimeStep(ts,0.0,.01);
392:   TSSetFromOptions(ts);

394:   direction[0] = direction[1] = 1;
395:   terminate[0] = terminate[1] = PETSC_FALSE;

397:   TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)ctx);

399:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
400:      Solve nonlinear system
401:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
402:   TSSolve(ts,U);

404:   TSGetSolveTime(ts,&ftime);
405:   TSGetTimeStepNumber(ts,&steps);

407:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
408:      Adjoint model starts here
409:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
410:   MatCreateVecs(A,&lambda[0],NULL);
411:   /*   Set initial conditions for the adjoint integration */
412:   VecGetArray(lambda[0],&y_ptr);
413:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
414:   VecRestoreArray(lambda[0],&y_ptr);

416:   MatCreateVecs(Jacp,&mu[0],NULL);
417:   VecGetArray(mu[0],&mx_ptr);
418:   mx_ptr[0] = -1.0;
419:   VecRestoreArray(mu[0],&mx_ptr);
420:   TSSetCostGradients(ts,1,lambda,mu);

422:   TSAdjointSetRHSJacobian(ts,Jacp,RHSJacobianP,ctx);

424:   TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
425:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
426:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_FALSE,ctx);

428:   TSAdjointSolve(ts);
429:   TSGetCostIntegral(ts,&q);
430:   ComputeSensiP(lambda[0],mu[0],ctx);

432:   VecGetArrayRead(q,&qx_ptr);
433:   *f   = -ctx->Pm + qx_ptr[0];
434:   VecRestoreArrayRead(q,&qx_ptr);

436:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
437:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
438:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
439:   MatDestroy(&A);
440:   MatDestroy(&Jacp);
441:   VecDestroy(&U);
442:   VecDestroy(&lambda[0]);
443:   VecDestroy(&mu[0]);
444:   TSDestroy(&ts);

446:   return 0;
447: }