Actual source code: ex3adj.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 <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: /*
 64:      Defines the ODE passed to the ODE solver
 65: */
 66: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
 67: {
 68:   PetscErrorCode    ierr;
 69:   PetscScalar       *f,Pmax;
 70:   const PetscScalar *u,*udot;

 73:   /*  The next three lines allow us to access the entries of the vectors directly */
 74:   VecGetArrayRead(U,&u);
 75:   VecGetArrayRead(Udot,&udot);
 76:   VecGetArray(F,&f);
 77:   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 */
 78:   else Pmax = ctx->Pmax;
 79: 
 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:   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 */
105:   else Pmax = ctx->Pmax;

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

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

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

125: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
126: {
128:   PetscInt       row[] = {0,1},col[]={0};
129:   PetscScalar    *x,J[2][1];
130: 
132:   VecGetArray(X,&x);

134:   J[0][0] = 0;
135:   J[1][0] = 1.;
136:   MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);

138:   MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
139:   MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
140:   return(0);
141: }

145: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
146: {
147:   PetscErrorCode    ierr;
148:   PetscScalar       *r;
149:   const PetscScalar *u;

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

162: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
163: {
164:   PetscErrorCode    ierr;
165:   PetscScalar       *ry;
166:   const PetscScalar *u;

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

179: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
180: {
181:   PetscErrorCode    ierr;
182:   PetscScalar       *rp;
183:   const PetscScalar *u;

186:   VecGetArrayRead(U,&u);
187:   VecGetArray(drdp[0],&rp);
188:   rp[0] = 0.;
189:   VecRestoreArray(drdp[0],&rp);
190:   VecGetArrayRead(U,&u);
191:   return(0);
192: }

196: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
197: {
198:   PetscErrorCode    ierr;
199:   PetscScalar       sensip;
200:   const PetscScalar *x,*y;
201: 
203:   VecGetArrayRead(lambda,&x);
204:   VecGetArrayRead(mu,&y);
205:   sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
206:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %.7f \n",(double)sensip);
207:   VecRestoreArrayRead(lambda,&x);
208:   VecRestoreArrayRead(mu,&y);
209:   return(0);
210: }

214: int main(int argc,char **argv)
215: {
216:   TS             ts;            /* ODE integrator */
217:   Vec            U;             /* solution will be stored here */
218:   Mat            A;             /* Jacobian matrix */
219:   Mat            Jacp;          /* Jacobian matrix */
221:   PetscMPIInt    size;
222:   PetscInt       n = 2;
223:   AppCtx         ctx;
224:   PetscScalar    *u;
225:   PetscReal      du[2] = {0.0,0.0};
226:   PetscBool      ensemble = PETSC_FALSE,flg1,flg2;
227:   PetscReal      ftime;
228:   PetscInt       steps;
229:   PetscScalar    *x_ptr,*y_ptr;
230:   Vec            lambda[1],q,mu[1];

232:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
233:      Initialize program
234:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
235:   PetscInitialize(&argc,&argv,(char*)0,help);
236:   MPI_Comm_size(PETSC_COMM_WORLD,&size);
237:   if (size > 1) SETERRQ(PETSC_COMM_WORLD,PETSC_ERR_SUP,"Only for sequential runs");

239:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
240:     Create necessary matrix and vectors
241:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
242:   MatCreate(PETSC_COMM_WORLD,&A);
243:   MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
244:   MatSetType(A,MATDENSE);
245:   MatSetFromOptions(A);
246:   MatSetUp(A);

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

250:   MatCreate(PETSC_COMM_WORLD,&Jacp);
251:   MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
252:   MatSetFromOptions(Jacp);
253:   MatSetUp(Jacp);

255:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
256:     Set runtime options
257:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
258:   PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
259:   {
260:     ctx.beta    = 2;
261:     ctx.c       = 10000.0;
262:     ctx.u_s     = 1.0;
263:     ctx.omega_s = 1.0;
264:     ctx.omega_b = 120.0*PETSC_PI;
265:     ctx.H       = 5.0;
266:     PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
267:     ctx.D       = 5.0;
268:     PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
269:     ctx.E       = 1.1378;
270:     ctx.V       = 1.0;
271:     ctx.X       = 0.545;
272:     ctx.Pmax    = ctx.E*ctx.V/ctx.X;;
273:     PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
274:     ctx.Pm      = 1.1;
275:     PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
276:     ctx.tf      = 0.1;
277:     ctx.tcl     = 0.2;
278:     PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
279:     PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
280:     PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);
281:     if (ensemble) {
282:       ctx.tf      = -1;
283:       ctx.tcl     = -1;
284:     }

286:     VecGetArray(U,&u);
287:     u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
288:     u[1] = 1.0;
289:     PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);
290:     n    = 2;
291:     PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);
292:     u[0] += du[0];
293:     u[1] += du[1];
294:     VecRestoreArray(U,&u);
295:     if (flg1 || flg2) {
296:       ctx.tf      = -1;
297:       ctx.tcl     = -1;
298:     }
299:   }
300:   PetscOptionsEnd();

302:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
303:      Create timestepping solver context
304:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
305:   TSCreate(PETSC_COMM_WORLD,&ts);
306:   TSSetProblemType(ts,TS_NONLINEAR);
307:   TSSetType(ts,TSBEULER);
308:   TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);
309:   TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);
310:   TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);

312:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
313:      Set initial conditions
314:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
315:   TSSetSolution(ts,U);

317:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
318:     Save trajectory of solution so that TSAdjointSolve() may be used
319:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
320:   TSSetSaveTrajectory(ts);

322:   MatCreateVecs(A,&lambda[0],NULL);
323:   /*   Set initial conditions for the adjoint integration */
324:   VecGetArray(lambda[0],&y_ptr);
325:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
326:   VecRestoreArray(lambda[0],&y_ptr);

328:   MatCreateVecs(Jacp,&mu[0],NULL);
329:   VecGetArray(mu[0],&x_ptr);
330:   x_ptr[0] = -1.0;
331:   VecRestoreArray(mu[0],&x_ptr);
332:   TSSetCostGradients(ts,1,lambda,mu);
333:   TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
334:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
335:                                         (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,&ctx);

337:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
338:      Set solver options
339:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
340:   TSSetDuration(ts,PETSC_DEFAULT,10.0);
341:   TSSetInitialTimeStep(ts,0.0,.01);
342:   TSSetFromOptions(ts);

344: #if 0
345:   TSSetPostStep(ts,PostStepFunction);
346: #endif
347:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
348:      Solve nonlinear system
349:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
350:   if (ensemble) {
351:     for (du[1] = -2.5; du[1] <= .01; du[1] += .1) {
352:       VecGetArray(U,&u);
353:       u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax);
354:       u[1] = ctx.omega_s;
355:       u[0] += du[0];
356:       u[1] += du[1];
357:       VecRestoreArray(U,&u);
358:       TSSetInitialTimeStep(ts,0.0,.01);
359:       TSSolve(ts,U);
360:     }
361:   } else {
362:     TSSolve(ts,U);
363:   }
364:   VecView(U,PETSC_VIEWER_STDOUT_WORLD);
365:   TSGetSolveTime(ts,&ftime);
366:   TSGetTimeStepNumber(ts,&steps);

368:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
369:      Adjoint model starts here
370:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
371:   /*   Set initial conditions for the adjoint integration */
372:   VecGetArray(lambda[0],&y_ptr);
373:   y_ptr[0] = 0.0; y_ptr[1] = 0.0;
374:   VecRestoreArray(lambda[0],&y_ptr);

376:   VecGetArray(mu[0],&x_ptr);
377:   x_ptr[0] = -1.0;
378:   VecRestoreArray(mu[0],&x_ptr);

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

383:   TSAdjointSolve(ts);

385:   PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt initial conditions: d[Psi(tf)]/d[phi0]  d[Psi(tf)]/d[omega0]\n");
386:   VecView(lambda[0],PETSC_VIEWER_STDOUT_WORLD);
387:   VecView(mu[0],PETSC_VIEWER_STDOUT_WORLD);
388:   TSGetCostIntegral(ts,&q);
389:   VecView(q,PETSC_VIEWER_STDOUT_WORLD);
390:   VecGetArray(q,&x_ptr);
391:   PetscPrintf(PETSC_COMM_WORLD,"\n cost function=%g\n",(double)(x_ptr[0]-ctx.Pm));
392: 
393:   ComputeSensiP(lambda[0],mu[0],&ctx);

395:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
396:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
397:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
398:   MatDestroy(&A);
399:   MatDestroy(&Jacp);
400:   VecDestroy(&U);
401:   VecDestroy(&lambda[0]);
402:   VecDestroy(&mu[0]);
403:   TSDestroy(&ts);
404:   PetscFinalize();
405:   return(0);
406: }