Actual source code: ex9opt.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: static char help[] = "Basic equation for generator stability analysis.\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}



Ensemble of initial conditions
./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly

Fault at .1 seconds
./ex2 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly

Initial conditions same as when fault is ended
./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly


 25: /*
 26:    Include "petscts.h" so that we can use TS solvers.  Note that this
 27:    file automatically includes:
 28:      petscsys.h       - base PETSc routines   petscvec.h - vectors
 29:      petscmat.h - matrices
 30:      petscis.h     - index sets            petscksp.h - Krylov subspace methods
 31:      petscviewer.h - viewers               petscpc.h  - preconditioners
 32:      petscksp.h   - linear solvers
 33: */
 34: #include <petsctao.h>
 35: #include <petscts.h>

 37: typedef struct {
 38:   PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c;
 39:   PetscInt    beta;
 40:   PetscReal   tf,tcl;
 41: } AppCtx;

 43: PetscErrorCode FormFunction(Tao,Vec,PetscReal*,void*);
 44: PetscErrorCode FormGradient(Tao,Vec,Vec,void*);

 48: /*
 49:      Defines the ODE passed to the ODE solver
 50: */
 51: static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx)
 52: {
 53:   PetscErrorCode    ierr;
 54:   PetscScalar       *f,Pmax;
 55:   const PetscScalar *u;

 58:   /*  The next three lines allow us to access the entries of the vectors directly */
 59:   VecGetArrayRead(U,&u);
 60:   VecGetArray(F,&f);
 61:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 62:   else Pmax = ctx->Pmax;

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

 67:   VecRestoreArrayRead(U,&u);
 68:   VecRestoreArray(F,&f);
 69:   return(0);
 70: }

 74: /*
 75:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
 76: */
 77: static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx)
 78: {
 79:   PetscErrorCode    ierr;
 80:   PetscInt          rowcol[] = {0,1};
 81:   PetscScalar       J[2][2],Pmax;
 82:   const PetscScalar *u;

 85:   VecGetArrayRead(U,&u);
 86:   if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */
 87:   else Pmax = ctx->Pmax;

 89:   J[0][0] = 0;                                  J[0][1] = ctx->omega_b;
 90:   J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H);  J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H);

 92:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
 93:   VecRestoreArrayRead(U,&u);

 95:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
 96:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
 97:   if (A != B) {
 98:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
 99:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
100:   }
101:   return(0);
102: }

106: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
107: {
109:   PetscInt       row[] = {0,1},col[]={0};
110:   PetscScalar    *x,J[2][1];
111:   AppCtx         *ctx=(AppCtx*)ctx0;

114:   VecGetArray(X,&x);

116:   J[0][0] = 0;
117:   J[1][0] = ctx->omega_s/(2.0*ctx->H);
118:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

120:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
121:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
122:   return(0);
123: }

127: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
128: {
129:   PetscErrorCode    ierr;
130:   PetscScalar       *r;
131:   const PetscScalar *u;

134:   VecGetArrayRead(U,&u);
135:   VecGetArray(R,&r);
136:   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
137:   VecRestoreArray(R,&r);
138:   VecRestoreArrayRead(U,&u);
139:   return(0);
140: }

144: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
145: {
146:   PetscErrorCode    ierr;
147:   PetscScalar       *ry;
148:   const PetscScalar *u;

151:   VecGetArrayRead(U,&u);
152:   VecGetArray(drdy[0],&ry);
153:   ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
154:   VecRestoreArray(drdy[0],&ry);
155:   VecRestoreArrayRead(U,&u);
156:   return(0);
157: }

161: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
162: {
164:   PetscScalar    *rp;

167:   VecGetArray(drdp[0],&rp);
168:   rp[0] = 0.;
169:   VecRestoreArray(drdp[0],&rp);
170:   return(0);
171: }

175: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
176: {
177:   PetscErrorCode    ierr;
178:   PetscScalar       *y,sensip;
179:   const PetscScalar *x;

182:   VecGetArrayRead(lambda,&x);
183:   VecGetArray(mu,&y);
184:   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
185:   /*PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %g \n",(double)sensip);*/
186:   y[0] = sensip;
187:   VecRestoreArray(mu,&y);
188:   VecRestoreArrayRead(lambda,&x);
189:   return(0);
190: }

194: int main(int argc,char **argv)
195: {
196:   Vec                p;
197:   PetscScalar        *x_ptr;
198:   PetscErrorCode     ierr;
199:   PetscMPIInt        size;
200:   AppCtx             ctx;
201:   Vec                lowerb,upperb;
202:   Tao                tao;
203:   TaoConvergedReason reason;
204:   KSP                ksp;
205:   PC                 pc;

207:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
208:      Initialize program
209:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
210:   PetscInitialize(&argc,&argv,NULL,help);
212:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
213:   if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");

215:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
216:     Set runtime options
217:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
218:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
219:   {
220:     ctx.beta    = 2;
221:     ctx.c       = 10000.0;
222:     ctx.u_s     = 1.0;
223:     ctx.omega_s = 1.0;
224:     ctx.omega_b = 120.0*PETSC_PI;
225:     ctx.H       = 5.0;
226:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
227:     ctx.D       = 5.0;
228:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
229:     ctx.E       = 1.1378;
230:     ctx.V       = 1.0;
231:     ctx.X       = 0.545;
232:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;;
233:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
234:     ctx.Pm      = 1.0194;
235:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
236:     ctx.tf      = 0.1;
237:     ctx.tcl     = 0.2;
238:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
239:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);

241:   }
242:   PetscOptionsEnd();

244:   /* Create TAO solver and set desired solution method */
245:   TaoCreate(PETSC_COMM_WORLD,&tao);
246:   TaoSetType(tao,TAOBLMVM);

248:   /*
249:      Optimization starts
250:   */
251:   /* Set initial solution guess */
252:   VecCreateSeq(PETSC_COMM_WORLD,1,&p);
253:   VecGetArray(p,&x_ptr);
254:   x_ptr[0]   = ctx.Pm;
255:   VecRestoreArray(p,&x_ptr);

257:   TaoSetInitialVector(tao,p);
258:   /* Set routine for function and gradient evaluation */
259:   TaoSetObjectiveRoutine(tao,FormFunction,(void *)&ctx);
260:   TaoSetGradientRoutine(tao,FormGradient,(void *)&ctx);

262:   /* Set bounds for the optimization */
263:   VecDuplicate(p,&lowerb);
264:   VecDuplicate(p,&upperb);
265:   VecGetArray(lowerb,&x_ptr);
266:   x_ptr[0] = 0.;
267:   VecRestoreArray(lowerb,&x_ptr);
268:   VecGetArray(upperb,&x_ptr);
269:   x_ptr[0] = 1.1;;
270:   VecRestoreArray(upperb,&x_ptr);
271:   TaoSetVariableBounds(tao,lowerb,upperb);

273:   /* Check for any TAO command line options */
274:   TaoSetFromOptions(tao);
275:   TaoGetKSP(tao,&ksp);
276:   if (ksp) {
277:     KSPGetPC(ksp,&pc);
278:     PCSetType(pc,PCNONE);
279:   }

281:   TaoSetTolerances(tao,1e-15,1e-15,1e-15,1e-15,1e-15);
282:   /* SOLVE THE APPLICATION */
283:   TaoSolve(tao);

285:   /* Get information on termination */
286:   TaoGetConvergedReason(tao,&reason);
287:   if (reason <= 0){
288:     ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
289:   }

291:   VecView(p,PETSC_VIEWER_STDOUT_WORLD);
292:   /* Free TAO data structures */
293:   TaoDestroy(&tao);
294:   VecDestroy(&p);
295:   VecDestroy(&lowerb);
296:   VecDestroy(&upperb);
297:   PetscFinalize();
298:   return 0;
299: }

301: /* ------------------------------------------------------------------ */
304: /*
305:    FormFunction - Evaluates the function

307:    Input Parameters:
308:    tao - the Tao context
309:    X   - the input vector
310:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()

312:    Output Parameters:
313:    f   - the newly evaluated function
314: */
315: PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
316: {
317:   AppCtx         *ctx = (AppCtx*)ctx0;
318:   TS             ts;

320:   Vec            U;             /* solution will be stored here */
321:   Mat            A;             /* Jacobian matrix */
323:   PetscInt       n = 2;
324:   PetscScalar    *u;
325:   PetscScalar    *x_ptr;
326:   Vec            lambda[1],q,mu[1];

328:   VecGetArray(P,&x_ptr);
329:   ctx->Pm = x_ptr[0];
330:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
331:     Create necessary matrix and vectors
332:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
333:   MatCreate(PETSC_COMM_WORLD,&A);
334:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
335:   MatSetType(A,MATDENSE);
336:   MatSetFromOptions(A);
337:   MatSetUp(A);

339:   MatCreateVecs(A,&U,NULL);
340:   MatCreateVecs(A,&lambda[0],NULL);
341:   MatCreateVecs(A,&mu[0],NULL);

343:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
344:      Create timestepping solver context
345:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
346:   TSCreate(PETSC_COMM_WORLD,&ts);
347:   TSSetProblemType(ts,TS_NONLINEAR);
348:   TSSetType(ts,TSRK);
349:   TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,ctx);
350:   TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,ctx);
351:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

353:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
354:      Set initial conditions
355:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
356:   VecGetArray(U,&u);
357:   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
358:   u[1] = 1.0;
359:   VecRestoreArray(U,&u);
360:   TSSetSolution(ts,U);

362:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
363:     Save trajectory of solution so that TSAdjointSolve() may be used
364:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
365:   TSSetSaveTrajectory(ts);

367:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
368:      Set solver options
369:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
370:   TSSetDuration(ts,PETSC_DEFAULT,10.0);
371:   TSSetInitialTimeStep(ts,0.0,.01);
372:   TSSetFromOptions(ts);

374:   TSSetCostGradients(ts,1,lambda,mu);
375:   TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
376:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
377:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,ctx);

379:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
380:      Solve nonlinear system
381:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
382:   TSSolve(ts,U);
383:   TSGetCostIntegral(ts,&q);
384:   VecGetArray(q,&x_ptr);
385:   *f   = -ctx->Pm + x_ptr[0];

387:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
388:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
389:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
390:   MatDestroy(&A);
391:   VecDestroy(&U);
392:   VecDestroy(&lambda[0]);
393:   VecDestroy(&mu[0]);
394:   TSDestroy(&ts);

396:   return 0;
397: }


400: PetscErrorCode FormGradient(Tao tao,Vec P,Vec G,void *ctx0)
401: {
402:   AppCtx         *ctx = (AppCtx*)ctx0;
403:   TS             ts;

405:   Vec            U;             /* solution will be stored here */
406:   Mat            A;             /* Jacobian matrix */
407:   Mat            Jacp;          /* Jacobian matrix */
409:   PetscInt       n = 2;
410:   PetscReal      ftime;
411:   PetscInt       steps;
412:   PetscScalar    *u;
413:   PetscScalar    *x_ptr,*y_ptr;
414:   Vec            lambda[1],q,mu[1];

416:   VecGetArray(P,&x_ptr);
417:   ctx->Pm = x_ptr[0];
418:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
419:     Create necessary matrix and vectors
420:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
421:   MatCreate(PETSC_COMM_WORLD,&A);
422:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
423:   MatSetType(A,MATDENSE);
424:   MatSetFromOptions(A);
425:   MatSetUp(A);

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

429:   MatCreate(PETSC_COMM_WORLD,&Jacp);
430:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
431:   MatSetFromOptions(Jacp);
432:   MatSetUp(Jacp);

434:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
435:      Create timestepping solver context
436:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
437:   TSCreate(PETSC_COMM_WORLD,&ts);
438:   TSSetProblemType(ts,TS_NONLINEAR);
439:   TSSetType(ts,TSRK);
440:   TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,ctx);
441:   TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,ctx);
442:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

444:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
445:      Set initial conditions
446:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
447:   VecGetArray(U,&u);
448:   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
449:   u[1] = 1.0;
450:   VecRestoreArray(U,&u);
451:   TSSetSolution(ts,U);

453:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
454:     Save trajectory of solution so that TSAdjointSolve() may be used
455:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
456:   TSSetSaveTrajectory(ts);

458:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
459:      Set solver options
460:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
461:   TSSetDuration(ts,PETSC_DEFAULT,10.0);
462:   TSSetInitialTimeStep(ts,0.0,.01);
463:   TSSetFromOptions(ts);

465:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
466:      Solve nonlinear system
467:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
468:   TSSolve(ts,U);

470:   TSGetSolveTime(ts,&ftime);
471:   TSGetTimeStepNumber(ts,&steps);

473:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
474:      Adjoint model starts here
475:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
476:   MatCreateVecs(A,&lambda[0],NULL);
477:   /*   Set initial conditions for the adjoint integration */
478:   VecGetArray(lambda[0],&y_ptr);
479:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
480:   VecRestoreArray(lambda[0],&y_ptr);

482:   MatCreateVecs(Jacp,&mu[0],NULL);
483:   VecGetArray(mu[0],&x_ptr);
484:   x_ptr[0] = -1.0;
485:   VecRestoreArray(mu[0],&x_ptr);
486:   TSSetCostGradients(ts,1,lambda,mu);

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

490:   TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
491:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
492:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,ctx);

494:   TSAdjointSolve(ts);
495:   TSGetCostIntegral(ts,&q);
496:   ComputeSensiP(lambda[0],mu[0],ctx);

498:   VecCopy(mu[0],G);

500:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
501:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
502:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
503:   MatDestroy(&A);
504:   MatDestroy(&Jacp);
505:   VecDestroy(&U);
506:   VecDestroy(&lambda[0]);
507:   VecDestroy(&mu[0]);
508:   TSDestroy(&ts);

510:   return 0;
511: }