Actual source code: ex9opt.c petsc-3.7.7 2017-09-25 
   
  2:  static char help[] = "Basic equation for generator stability analysis.\n" ;
 
                 
 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:  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:  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:  PetscErrorCode  RHSJacobianP(TS  ts,PetscReal  t,Vec  X,Mat  A,void *ctx0)107:  {
109:    PetscInt        row[] = {0,1},col[]={0};
110:    PetscScalar     J[2][1];
111:    AppCtx         *ctx=(AppCtx*)ctx0;
114:    J[0][0] = 0;
115:    J[1][0] = ctx->omega_s/(2.0*ctx->H);
116:    MatSetValues (A,2,row,1,col,&J[0][0],INSERT_VALUES );
117:    MatAssemblyBegin (A,MAT_FINAL_ASSEMBLY);
118:    MatAssemblyEnd (A,MAT_FINAL_ASSEMBLY);
119:    return (0);
120:  }
124:  PetscErrorCode  CostIntegrand(TS  ts,PetscReal  t,Vec  U,Vec  R,AppCtx *ctx)125:  {
126:    PetscErrorCode     ierr;
127:    PetscScalar        *r;
128:    const PetscScalar  *u;
131:    VecGetArrayRead (U,&u);
132:    VecGetArray (R,&r);
133:    r[0] = ctx->c*PetscPowScalarInt(PetscMax (0., u[0]-ctx->u_s),ctx->beta);
134:    VecRestoreArray (R,&r);
135:    VecRestoreArrayRead (U,&u);
136:    return (0);
137:  }
141:  PetscErrorCode  DRDYFunction(TS  ts,PetscReal  t,Vec  U,Vec  *drdy,AppCtx *ctx)142:  {
143:    PetscErrorCode     ierr;
144:    PetscScalar        *ry;
145:    const PetscScalar  *u;
148:    VecGetArrayRead (U,&u);
149:    VecGetArray (drdy[0],&ry);
150:    ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax (0., u[0]-ctx->u_s),ctx->beta-1);
151:    VecRestoreArray (drdy[0],&ry);
152:    VecRestoreArrayRead (U,&u);
153:    return (0);
154:  }
158:  PetscErrorCode  DRDPFunction(TS  ts,PetscReal  t,Vec  U,Vec  *drdp,AppCtx *ctx)159:  {
161:    PetscScalar     *rp;
164:    VecGetArray (drdp[0],&rp);
165:    rp[0] = 0.;
166:    VecRestoreArray (drdp[0],&rp);
167:    return (0);
168:  }
172:  PetscErrorCode  ComputeSensiP(Vec  lambda,Vec  mu,AppCtx *ctx)173:  {
174:    PetscErrorCode     ierr;
175:    PetscScalar        *y,sensip;
176:    const PetscScalar  *x;
179:    VecGetArrayRead (lambda,&x);
180:    VecGetArray (mu,&y);
181:    sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
182:    /*PetscPrintf (PETSC_COMM_WORLD ,"\n sensitivity wrt parameter pm: %g \n",(double)sensip);*/ 
183:    y[0] = sensip;
184:    VecRestoreArray (mu,&y);
185:    VecRestoreArrayRead (lambda,&x);
186:    return (0);
187:  }
191:  192:  {
193:    Vec                 p;
194:    PetscScalar         *x_ptr;
195:    PetscErrorCode      ierr;
196:    PetscMPIInt         size;
197:    AppCtx             ctx;
198:    Vec                 lowerb,upperb;
199:    Tao                 tao;
200:    KSP                 ksp;
201:    PC                  pc;
203:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
204:       Initialize program 
205:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
206:    PetscInitialize (&argc,&argv,NULL,help);
208:    MPI_Comm_size (PETSC_COMM_WORLD ,&size);
209:    if  (size != 1) SETERRQ (PETSC_COMM_SELF ,1,"This is a uniprocessor example only!" );
211:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
212:      Set runtime options 
213:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
214:    PetscOptionsBegin (PETSC_COMM_WORLD ,NULL,"Swing equation options" ,"" );
215:    {
216:      ctx.beta    = 2;
217:      ctx.c       = 10000.0;
218:      ctx.u_s     = 1.0;
219:      ctx.omega_s = 1.0;
220:      ctx.omega_b = 120.0*PETSC_PI;
221:      ctx.H       = 5.0;
222:      PetscOptionsScalar ("-Inertia" ,"" ,"" ,ctx.H,&ctx.H,NULL);
223:      ctx.D       = 5.0;
224:      PetscOptionsScalar ("-D" ,"" ,"" ,ctx.D,&ctx.D,NULL);
225:  #if defined(PETSC_USE_REAL___FLOAT128) 
226:      ctx.E       = 1.1378q;
227:  #else 
228:      ctx.E       = 1.1378;
229:  #endif 
230:      ctx.V       = 1.0;
231:  #if defined(PETSC_USE_REAL___FLOAT128) 
232:      ctx.X       = 0.545q;
233:  #else 
234:      ctx.X       = 0.545;
235:  #endif 
236:      ctx.Pmax    = ctx.E*ctx.V/ctx.X;;
237:      PetscOptionsScalar ("-Pmax" ,"" ,"" ,ctx.Pmax,&ctx.Pmax,NULL);
238:  #if defined(PETSC_USE_REAL___FLOAT128) 
239:      ctx.Pm      = 1.0194q;
240:  #else 
241:      ctx.Pm      = 1.0194;
242:  #endif 
243:      PetscOptionsScalar ("-Pm" ,"" ,"" ,ctx.Pm,&ctx.Pm,NULL);
244:  #if defined(PETSC_USE_REAL___FLOAT128) 
245:      ctx.tf      = 0.1q;
246:      ctx.tcl     = 0.2q;
247:  #else 
248:      ctx.tf      = 0.1;
249:      ctx.tcl     = 0.2;
250:  #endif 
251:      PetscOptionsReal ("-tf" ,"Time to start fault" ,"" ,ctx.tf,&ctx.tf,NULL);
252:      PetscOptionsReal ("-tcl" ,"Time to end fault" ,"" ,ctx.tcl,&ctx.tcl,NULL);
254:    }
255:    PetscOptionsEnd ();
257:    /* Create TAO solver and set desired solution method */ 
258:    TaoCreate (PETSC_COMM_WORLD ,&tao);
259:    TaoSetType (tao,TAOBLMVM );
261:    /* 
262:       Optimization starts 
263:    */ 
264:    /* Set initial solution guess */ 
265:    VecCreateSeq (PETSC_COMM_WORLD ,1,&p);
266:    VecGetArray (p,&x_ptr);
267:    x_ptr[0]   = ctx.Pm;
268:    VecRestoreArray (p,&x_ptr);
270:    TaoSetInitialVector (tao,p);
271:    /* Set routine for function and gradient evaluation */ 
272:    TaoSetObjectiveRoutine (tao,FormFunction,(void *)&ctx);
273:    TaoSetGradientRoutine (tao,FormGradient,(void *)&ctx);
275:    /* Set bounds for the optimization */ 
276:    VecDuplicate (p,&lowerb);
277:    VecDuplicate (p,&upperb);
278:    VecGetArray (lowerb,&x_ptr);
279:    x_ptr[0] = 0.;
280:    VecRestoreArray (lowerb,&x_ptr);
281:    VecGetArray (upperb,&x_ptr);
282:  #if defined(PETSC_USE_REAL___FLOAT128) 
283:    x_ptr[0] = 1.1q;
284:  #else 
285:    x_ptr[0] = 1.1;
286:  #endif 
287:    VecRestoreArray (upperb,&x_ptr);
288:    TaoSetVariableBounds (tao,lowerb,upperb);
290:    /* Check for any TAO command line options */ 
291:    TaoSetFromOptions (tao);
292:    TaoGetKSP (tao,&ksp);
293:    if  (ksp) {
294:      KSPGetPC (ksp,&pc);
295:      PCSetType (pc,PCNONE );
296:    }
298:    /* SOLVE THE APPLICATION */ 
299:    TaoSolve (tao);
301:    VecView (p,PETSC_VIEWER_STDOUT_WORLD );
302:    /* Free TAO data structures */ 
303:    TaoDestroy (&tao);
304:    VecDestroy (&p);
305:    VecDestroy (&lowerb);
306:    VecDestroy (&upperb);
307:    PetscFinalize ();
308:    return  0;
309:  }
311:  /* ------------------------------------------------------------------ */ 
314:  /* 
315:     FormFunction - Evaluates the function 
317:     Input Parameters: 
318:     tao - the Tao  context 
319:     X   - the input vector 
320:     ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine () 
322:     Output Parameters: 
323:     f   - the newly evaluated function 
324:  */ 
325:  PetscErrorCode  FormFunction(Tao  tao,Vec  P,PetscReal  *f,void *ctx0)326:  {
327:    AppCtx         *ctx = (AppCtx*)ctx0;
328:    TS              ts;
330:    Vec             U;             /* solution will be stored here */ 
331:    Mat             A;             /* Jacobian matrix */ 
333:    PetscInt        n = 2;
334:    PetscScalar     *u;
335:    PetscScalar     *x_ptr;
336:    Vec             lambda[1],q,mu[1];
338:    VecGetArray (P,&x_ptr);
339:    ctx->Pm = x_ptr[0];
340:    VecRestoreArray (P,&x_ptr);
341:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
342:      Create necessary matrix and vectors 
343:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
344:    MatCreate (PETSC_COMM_WORLD ,&A);
345:    MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
346:    MatSetType (A,MATDENSE );
347:    MatSetFromOptions (A);
348:    MatSetUp (A);
350:    MatCreateVecs (A,&U,NULL);
351:    MatCreateVecs (A,&lambda[0],NULL);
352:    MatCreateVecs (A,&mu[0],NULL);
354:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
355:       Create timestepping solver context 
356:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
357:    TSCreate (PETSC_COMM_WORLD ,&ts);
358:    TSSetProblemType (ts,TS_NONLINEAR);
359:    TSSetType (ts,TSRK );
360:    TSSetRHSFunction (ts,NULL,(TSRHSFunction)RHSFunction,ctx);
361:    TSSetRHSJacobian (ts,A,A,(TSRHSJacobian)RHSJacobian,ctx);
362:    TSSetExactFinalTime (ts,TS_EXACTFINALTIME_MATCHSTEP);
364:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
365:       Set initial conditions 
366:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
367:    VecGetArray (U,&u);
368:    u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
369:    u[1] = 1.0;
370:    VecRestoreArray (U,&u);
371:    TSSetSolution (ts,U);
373:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
374:      Save trajectory of solution so that TSAdjointSolve () may be used 
375:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
376:    TSSetSaveTrajectory (ts);
378:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
379:       Set solver options 
380:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
381:    TSSetDuration (ts,PETSC_DEFAULT ,1.0);
382:  #if defined(PETSC_USE_REAL___FLOAT128) 
383:    TSSetInitialTimeStep (ts,0.0,.01q);
384:  #else 
385:    TSSetInitialTimeStep (ts,0.0,.01);
386:  #endif 
387:    TSSetFromOptions (ts);
389:    TSSetCostGradients (ts,1,lambda,mu);
390:    TSSetCostIntegrand (ts,1,(PetscErrorCode  (*)(TS ,PetscReal ,Vec ,Vec ,void*))CostIntegrand,
391:                                          (PetscErrorCode  (*)(TS ,PetscReal ,Vec ,Vec *,void*))DRDYFunction,
392:                                          (PetscErrorCode  (*)(TS ,PetscReal ,Vec ,Vec *,void*))DRDPFunction,PETSC_TRUE ,ctx);
394:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
395:       Solve nonlinear system 
396:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
397:    TSSolve (ts,U);
398:    TSGetCostIntegral (ts,&q);
399:    VecGetArray (q,&x_ptr);
400:    *f   = -ctx->Pm + x_ptr[0];
401:    VecRestoreArray (q,&x_ptr);
403:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
404:       Free work space.  All PETSc objects should be destroyed when they are no longer needed. 
405:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
406:    MatDestroy (&A);
407:    VecDestroy (&U);
408:    VecDestroy (&lambda[0]);
409:    VecDestroy (&mu[0]);
410:    TSDestroy (&ts);
412:    return  0;
413:  }
416:  PetscErrorCode  FormGradient(Tao  tao,Vec  P,Vec  G,void *ctx0)417:  {
418:    AppCtx         *ctx = (AppCtx*)ctx0;
419:    TS              ts;
421:    Vec             U;             /* solution will be stored here */ 
422:    Mat             A;             /* Jacobian matrix */ 
423:    Mat             Jacp;          /* Jacobian matrix */ 
425:    PetscInt        n = 2;
426:    PetscReal       ftime;
427:    PetscInt        steps;
428:    PetscScalar     *u;
429:    PetscScalar     *x_ptr,*y_ptr;
430:    Vec             lambda[1],q,mu[1];
432:    VecGetArray (P,&x_ptr);
433:    ctx->Pm = x_ptr[0];
434:    VecRestoreArray (P,&x_ptr);
436:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
437:      Create necessary matrix and vectors 
438:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
439:    MatCreate (PETSC_COMM_WORLD ,&A);
440:    MatSetSizes (A,n,n,PETSC_DETERMINE ,PETSC_DETERMINE );
441:    MatSetType (A,MATDENSE );
442:    MatSetFromOptions (A);
443:    MatSetUp (A);
445:    MatCreateVecs (A,&U,NULL);
447:    MatCreate (PETSC_COMM_WORLD ,&Jacp);
448:    MatSetSizes (Jacp,PETSC_DECIDE ,PETSC_DECIDE ,2,1);
449:    MatSetFromOptions (Jacp);
450:    MatSetUp (Jacp);
452:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
453:       Create timestepping solver context 
454:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
455:    TSCreate (PETSC_COMM_WORLD ,&ts);
456:    TSSetProblemType (ts,TS_NONLINEAR);
457:    TSSetType (ts,TSRK );
458:    TSSetRHSFunction (ts,NULL,(TSRHSFunction)RHSFunction,ctx);
459:    TSSetRHSJacobian (ts,A,A,(TSRHSJacobian)RHSJacobian,ctx);
460:    TSSetExactFinalTime (ts,TS_EXACTFINALTIME_MATCHSTEP);
462:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
463:       Set initial conditions 
464:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
465:    VecGetArray (U,&u);
466:    u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
467:    u[1] = 1.0;
468:    VecRestoreArray (U,&u);
469:    TSSetSolution (ts,U);
471:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
472:      Save trajectory of solution so that TSAdjointSolve () may be used 
473:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
474:    TSSetSaveTrajectory (ts);
476:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
477:       Set solver options 
478:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
479:    TSSetDuration (ts,PETSC_DEFAULT ,1.0);
480:  #if defined(PETSC_USE_REAL___FLOAT128) 
481:    TSSetInitialTimeStep (ts,0.0,.01q);
482:  #else 
483:    TSSetInitialTimeStep (ts,0.0,.01);
484:  #endif 
485:    TSSetFromOptions (ts);
487:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
488:       Solve nonlinear system 
489:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
490:    TSSolve (ts,U);
492:    TSGetSolveTime (ts,&ftime);
493:    TSGetTimeStepNumber (ts,&steps);
495:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
496:       Adjoint model starts here 
497:       - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
498:    MatCreateVecs (A,&lambda[0],NULL);
499:    /*   Set initial conditions for the adjoint integration */ 
500:    VecGetArray (lambda[0],&y_ptr);
501:    y_ptr[0] = 0.0; y_ptr[1] = 0.0;
502:    VecRestoreArray (lambda[0],&y_ptr);
504:    MatCreateVecs (Jacp,&mu[0],NULL);
505:    VecGetArray (mu[0],&x_ptr);
506:    x_ptr[0] = -1.0;
507:    VecRestoreArray (mu[0],&x_ptr);
508:    TSSetCostGradients (ts,1,lambda,mu);
510:    TSAdjointSetRHSJacobian (ts,Jacp,RHSJacobianP,ctx);
512:    TSSetCostIntegrand (ts,1,(PetscErrorCode  (*)(TS ,PetscReal ,Vec ,Vec ,void*))CostIntegrand,
513:                                          (PetscErrorCode  (*)(TS ,PetscReal ,Vec ,Vec *,void*))DRDYFunction,
514:                                          (PetscErrorCode  (*)(TS ,PetscReal ,Vec ,Vec *,void*))DRDPFunction,PETSC_FALSE ,ctx);
516:    TSAdjointSolve (ts);
517:    TSGetCostIntegral (ts,&q);
518:    ComputeSensiP(lambda[0],mu[0],ctx);
520:    VecCopy (mu[0],G);
522:    /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 
523:       Free work space.  All PETSc objects should be destroyed when they are no longer needed. 
524:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 
525:    MatDestroy (&A);
526:    MatDestroy (&Jacp);
527:    VecDestroy (&U);
528:    VecDestroy (&lambda[0]);
529:    VecDestroy (&mu[0]);
530:    TSDestroy (&ts);
532:    return  0;
533:  }