Actual source code: ex3adj_events.c

petsc-3.6.4 2016-04-12
Report Typos and Errors
  2: static char help[] = "This is a version of ex3adj.c that also uses event monitor (TSEvent).\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 <petscts.h>

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

 44: PetscErrorCode PostStepFunction(TS ts)
 45: {
 47:   Vec            U;
 48:   PetscReal      t;
 49:   const PetscScalar    *u;

 52:   TSGetTime(ts,&t);
 53:   TSGetSolution(ts,&U);
 54:   VecGetArrayRead(U,&u);
 55:   PetscPrintf(PETSC_COMM_SELF,"delta(%3.2f) = %8.7f\n",t,u[0]);
 56:   VecRestoreArrayRead(U,&u);
 57: 
 58:   return(0);
 59: }

 63: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec U,PetscScalar *fvalue,void *app)
 64: {
 65:   AppCtx            *ctx=(AppCtx*)app;

 68:   /* Event for fault-on time */
 69:   fvalue[0] = t - ctx->tf;
 70:   /* Event for fault-off time */
 71:   fvalue[1] = t - ctx->tcl;
 72:   return(0);
 73: }

 77: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec U,PetscBool forwardsolve,void* app)
 78: {
 80:   AppCtx         *ctx=(AppCtx*)app;

 83:   if(forwardsolve) {
 84:     if (event_list[0] == 0) { /* Fault-on */
 85:       PetscPrintf(PETSC_COMM_SELF,"Applying fault at time = %3.2f\n",t);
 86:       /* Set Pmax to 0 */
 87:       ctx->Pmax = 0.0;
 88:     } else if (event_list[0] == 1) { /* Fault-off */
 89:       /* Set Pmax to initial value */
 90:       PetscPrintf(PETSC_COMM_SELF,"Removing fault at time = %3.2f\n",t);
 91:       ctx->Pmax = ctx->E*ctx->V/ctx->X;
 92:     }
 93:   } else { /* Reverse the events for the adjoint solve */
 94:     if (event_list[0] == 1) { /* Fault-on */
 95:       PetscPrintf(PETSC_COMM_SELF,"Adjoint solve:Applying fault at time = %3.2f\n",t);
 96:       /* Set Pmax to 0 */
 97:       ctx->Pmax = 0.0;
 98:     } else if (event_list[0] == 0) { /* Fault-off */
 99:       /* Set Pmax to initial value */
100:       PetscPrintf(PETSC_COMM_SELF,"Adjoint solve:Removing fault at time = %3.2f\n",t);
101:       ctx->Pmax = ctx->E*ctx->V/ctx->X;
102:     }
103:   }
104: 
105:   return(0);
106: }

110: /*
111:      Defines the ODE passed to the ODE solver
112: */
113: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
114: {
115:   PetscErrorCode    ierr;
116:   PetscScalar       *f;
117:   const PetscScalar *u,*udot;

120:   /*  The next three lines allow us to access the entries of the vectors directly */
121:   VecGetArrayRead(U,&u);
122:   VecGetArrayRead(Udot,&udot);
123:   VecGetArray(F,&f);
124: 
125:   f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s);
126:   f[1] = 2.0*ctx->H/ctx->omega_s*udot[1] +  ctx->Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - ctx->omega_s)- ctx->Pm;

128:   VecRestoreArrayRead(U,&u);
129:   VecRestoreArrayRead(Udot,&udot);
130:   VecRestoreArray(F,&f);
131:   return(0);
132: }

136: /*
137:      Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
138: */
139: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
140: {
141:   PetscErrorCode    ierr;
142:   PetscInt          rowcol[] = {0,1};
143:   PetscScalar       J[2][2];
144:   const PetscScalar *u,*udot;

147:   VecGetArrayRead(U,&u);
148:   VecGetArrayRead(Udot,&udot);

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

153:   MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
154:   VecRestoreArrayRead(U,&u);
155:   VecRestoreArrayRead(Udot,&udot);

157:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
158:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
159:   if (A != B) {
160:     MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
161:     MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
162:   }
163:   return(0);
164: }

168: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
169: {
171:   PetscInt       row[] = {0,1},col[]={0};
172:   PetscScalar    *x,J[2][1];
173: 
175:   VecGetArray(X,&x);

177:   J[0][0] = 0;
178:   J[1][0] = 1.;
179:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

181:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
182:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
183:   return(0);
184: }

188: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
189: {
190:   PetscErrorCode    ierr;
191:   PetscScalar       *r;
192:   const PetscScalar *u;

195:   VecGetArrayRead(U,&u);
196:   VecGetArray(R,&r);
197:   r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
198:   VecRestoreArray(R,&r);
199:   VecRestoreArrayRead(U,&u);
200:   return(0);
201: }

205: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
206: {
207:   PetscErrorCode    ierr;
208:   PetscScalar       *ry;
209:   const PetscScalar *u;

212:   VecGetArrayRead(U,&u);
213:   VecGetArray(drdy[0],&ry);
214:   ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
215:   VecRestoreArray(drdy[0],&ry);
216:   VecRestoreArrayRead(U,&u);
217:   return(0);
218: }

222: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
223: {
224:   PetscErrorCode    ierr;
225:   PetscScalar       *rp;
226:   const PetscScalar *u;

229:   VecGetArrayRead(U,&u);
230:   VecGetArray(drdp[0],&rp);
231:   rp[0] = 0.;
232:   VecRestoreArray(drdp[0],&rp);
233:   VecGetArrayRead(U,&u);
234:   return(0);
235: }

239: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
240: {
241:   PetscErrorCode    ierr;
242:   PetscScalar       sensip;
243:   const PetscScalar *x,*y;
244: 
246:   VecGetArrayRead(lambda,&x);
247:   VecGetArrayRead(mu,&y);
248:   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
249:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);
250:   VecRestoreArrayRead(lambda,&x);
251:   VecRestoreArrayRead(mu,&y);
252:   return(0);
253: }

257: int main(int argc,char **argv)
258: {
259:   TS             ts;            /* ODE integrator */
260:   Vec            U;             /* solution will be stored here */
261:   Mat            A;             /* Jacobian matrix */
262:   Mat            Jacp;          /* Jacobian matrix */
264:   PetscMPIInt    size;
265:   PetscInt       n = 2;
266:   AppCtx         ctx;
267:   PetscScalar    *u;
268:   PetscReal      du[2] = {0.0,0.0};
269:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
270:   PetscReal      ftime;
271:   PetscInt       steps;
272:   PetscScalar    *x_ptr,*y_ptr;
273:   Vec            lambda[1],q,mu[1];
274:   PetscInt       direction[2];
275:   PetscBool      terminate[2];


278:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
279:      Initialize program
280:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
281:   PetscInitialize(&argc,&argv,(char*)0,help);
282:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
283:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

285:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
286:     Create necessary matrix and vectors
287:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
288:   MatCreate(PETSC_COMM_WORLD,&A);
289:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
290:   MatSetType(A,MATDENSE);
291:   MatSetFromOptions(A);
292:   MatSetUp(A);

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

296:   MatCreate(PETSC_COMM_WORLD,&Jacp);
297:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
298:   MatSetFromOptions(Jacp);
299:   MatSetUp(Jacp);

301:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
302:     Set runtime options
303:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
304:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
305:   {
306:     ctx.beta    = 2;
307:     ctx.c       = 10000.0;
308:     ctx.u_s     = 1.0;
309:     ctx.omega_s = 1.0;
310:     ctx.omega_b = 120.0*PETSC_PI;
311:     ctx.H       = 5.0;
312:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
313:     ctx.D       = 5.0;
314:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
315:     ctx.E       = 1.1378;
316:     ctx.V       = 1.0;
317:     ctx.X       = 0.545;
318:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;;
319:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
320:     ctx.Pm      = 1.1;
321:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
322:     ctx.tf      = 0.1;
323:     ctx.tcl     = 0.2;
324:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
325:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
326:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
327:     if (ensemble) {
328:       ctx.tf      = -1;
329:       ctx.tcl     = -1;
330:     }

332:     VecGetArray(U,&u);
333:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
334:     u[1] = 1.0;
335:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
336:     n    = 2;
337:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
338:     u[0] += du[0];
339:     u[1] += du[1];
340:     VecRestoreArray(U,&u);
341:     if (flg1 || flg2) {
342:       ctx.tf      = -1;
343:       ctx.tcl     = -1;
344:     }
345:   }
346:   PetscOptionsEnd();

348:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
349:      Create timestepping solver context
350:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
351:   TSCreate(PETSC_COMM_WORLD,&ts);
352:   TSSetProblemType(ts,TS_NONLINEAR);
353:   TSSetType(ts,TSBEULER);
354:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
355:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);
356:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

358:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
359:      Set initial conditions
360:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
361:   TSSetSolution(ts,U);

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

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

375: #if 0
376:   TSSetPostStep(ts,PostStepFunction);
377: #endif
378:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
379:      Set events
380:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
381:   direction[0] = direction[1] = 1;
382:   terminate[0] = terminate[1] = PETSC_FALSE;
383: 
384:   TSSetEventMonitor(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)&ctx);

386:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
387:      Solve nonlinear system
388:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
389:   if (ensemble) {
390:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
391:       VecGetArray(U,&u);
392:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
393:       u[1] = ctx.omega_s;
394:       u[0] += du[0];
395:       u[1] += du[1];
396:       VecRestoreArray(U,&u);
397:       TSSetInitialTimeStep(ts,0.0,.01);
398:       TSSolve(ts,U);
399:     }
400:   } else {
401:     TSSolve(ts,U);
402:   }

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],&x_ptr);
418:   x_ptr[0] = -1.0;
419:   VecRestoreArray(mu[0],&x_ptr);
420:   TSSetCostGradients(ts,1,lambda,mu);

422:   /*   Set RHS JacobianP */
423:   TSAdjointSetRHSJacobian(ts,Jacp,RHSJacobianP,&ctx);

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

429:   TSAdjointSolve(ts);

431:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");
432:   VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
433:   VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
434:   TSGetCostIntegral(ts,&q);
435:   VecView(q,PETSC_VIEWER_STDOUT_WORLD);
436:   VecGetArray(q,&x_ptr);
437:   PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));
438: 
439:   ComputeSensiP(lambda[0],mu[0],&ctx);

441:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
442:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
443:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
444:   MatDestroy(&A);
445:   MatDestroy(&Jacp);
446:   VecDestroy(&U);
447:   VecDestroy(&lambda[0]);
448:   VecDestroy(&mu[0]);
449:   TSDestroy(&ts);
450:   PetscFinalize();
451:   return(0);
452: }