Actual source code: ex3opt.c
petsc-3.6.4 2016-04-12
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}
13: /*
14: Include "petscts.h" so that we can use TS solvers. Note that this
15: file automatically includes:
16: petscsys.h - base PETSc routines petscvec.h - vectors
17: petscmat.h - matrices
18: petscis.h - index sets petscksp.h - Krylov subspace methods
19: petscviewer.h - viewers petscpc.h - preconditioners
20: petscksp.h - linear solvers
21: */
22: #include <petsctao.h>
23: #include <petscts.h>
25: typedef struct {
26: PetscScalar H,D,omega_b,omega_s,Pmax,Pm,E,V,X,u_s,c;
27: PetscInt beta;
28: PetscReal tf,tcl;
29: } AppCtx;
31: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
35: /*
36: Defines the ODE passed to the ODE solver
37: */
38: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
39: {
40: PetscErrorCode ierr;
41: PetscScalar *f,Pmax;
42: const PetscScalar *u,*udot;
45: /* The next three lines allow us to access the entries of the vectors directly */
46: VecGetArrayRead(U,&u);
47: VecGetArrayRead(Udot,&udot);
48: VecGetArray(F,&f);
49: 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 */
50: else Pmax = ctx->Pmax;
52: f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s);
53: f[1] = 2.0*ctx->H/ctx->omega_s*udot[1] + Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - ctx->omega_s)- ctx->Pm;
55: VecRestoreArrayRead(U,&u);
56: VecRestoreArrayRead(Udot,&udot);
57: VecRestoreArray(F,&f);
58: return(0);
59: }
63: /*
64: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
65: */
66: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
67: {
68: PetscErrorCode ierr;
69: PetscInt rowcol[] = {0,1};
70: PetscScalar J[2][2],Pmax;
71: const PetscScalar *u,*udot;
74: VecGetArrayRead(U,&u);
75: VecGetArrayRead(Udot,&udot);
76: 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 */
77: else Pmax = ctx->Pmax;
79: J[0][0] = a; J[0][1] = -ctx->omega_b;
80: J[1][1] = 2.0*ctx->H/ctx->omega_s*a + ctx->D; J[1][0] = Pmax*PetscCosScalar(u[0]);
82: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
83: VecRestoreArrayRead(U,&u);
84: VecRestoreArrayRead(Udot,&udot);
86: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
87: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
88: if (A != B) {
89: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
90: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
91: }
92: return(0);
93: }
97: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
98: {
100: PetscInt row[] = {0,1},col[]={0};
101: PetscScalar J[2][1];
104: J[0][0] = 0;
105: J[1][0] = 1.;
106: MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
107: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
108: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
109: return(0);
110: }
114: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
115: {
116: PetscErrorCode ierr;
117: PetscScalar *r;
118: const PetscScalar *u;
121: VecGetArrayRead(U,&u);
122: VecGetArray(R,&r);
123: r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
124: VecRestoreArray(R,&r);
125: VecRestoreArrayRead(U,&u);
126: return(0);
127: }
131: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
132: {
133: PetscErrorCode ierr;
134: PetscScalar *ry;
135: const PetscScalar *u;
138: VecGetArrayRead(U,&u);
139: VecGetArray(drdy[0],&ry);
140: ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
141: VecRestoreArray(drdy[0],&ry);
142: VecRestoreArrayRead(U,&u);
143: return(0);
144: }
148: static PetscErrorCode DRDPFunction(TS ts,PetscReal t,Vec U,Vec *drdp,AppCtx *ctx)
149: {
151: PetscScalar *rp;
154: VecGetArray(drdp[0],&rp);
155: rp[0] = 0.;
156: VecRestoreArray(drdp[0],&rp);
157: return(0);
158: }
162: PetscErrorCode PostStep(TS ts)
163: {
165: Vec X;
166: PetscReal t;
169: TSGetTime(ts,&t);
170: if (t >= .2) {
171: TSGetSolution(ts,&X);
172: /* VecView(X,PETSC_VIEWER_STDOUT_WORLD);*/
173: exit(0);
174: /* results in initial conditions after fault of -u 0.496792,1.00932 */
175: }
176: return(0);
177: }
181: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
182: {
183: PetscErrorCode ierr;
184: PetscScalar *y,sensip;
185: const PetscScalar *x;
188: VecGetArrayRead(lambda,&x);
189: VecGetArray(mu,&y);
190: sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
191: /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %g \n",(double)sensip); */
192: y[0] = sensip;
193: VecRestoreArray(mu,&y);
194: VecRestoreArrayRead(lambda,&x);
195: return(0);
196: }
200: int main(int argc,char **argv)
201: {
202: Vec p;
203: PetscScalar *x_ptr;
204: PetscErrorCode ierr;
205: PetscMPIInt size;
206: AppCtx ctx;
207: Tao tao;
208: TaoConvergedReason reason;
209: KSP ksp;
210: PC pc;
211: Vec lowerb,upperb;
213: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
214: Initialize program
215: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
216: PetscInitialize(&argc,&argv,NULL,help);
218: MPI_Comm_size(PETSC_COMM_WORLD,&size);
219: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
221: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
222: Set runtime options
223: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
224: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
225: {
226: ctx.beta = 2;
227: ctx.c = 10000.0;
228: ctx.u_s = 1.0;
229: ctx.omega_s = 1.0;
230: ctx.omega_b = 120.0*PETSC_PI;
231: ctx.H = 5.0;
232: PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
233: ctx.D = 5.0;
234: PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
235: ctx.E = 1.1378;
236: ctx.V = 1.0;
237: ctx.X = 0.545;
238: ctx.Pmax = ctx.E*ctx.V/ctx.X;;
239: PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
240: ctx.Pm = 1.06161;
241: PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
242: ctx.tf = 0.1;
243: ctx.tcl = 0.2;
244: PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
245: PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
247: }
248: PetscOptionsEnd();
250: /* Create TAO solver and set desired solution method */
251: TaoCreate(PETSC_COMM_WORLD,&tao);
252: TaoSetType(tao,TAOBLMVM);
254: /*
255: Optimization starts
256: */
257: /* Set initial solution guess */
258: VecCreateSeq(PETSC_COMM_WORLD,1,&p);
259: VecGetArray(p,&x_ptr);
260: x_ptr[0] = ctx.Pm;
261: VecRestoreArray(p,&x_ptr);
263: TaoSetInitialVector(tao,p);
264: /* Set routine for function and gradient evaluation */
265: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(void *)&ctx);
267: /* Set bounds for the optimization */
268: VecDuplicate(p,&lowerb);
269: VecDuplicate(p,&upperb);
270: VecGetArray(lowerb,&x_ptr);
271: x_ptr[0] = 0.;
272: VecRestoreArray(lowerb,&x_ptr);
273: VecGetArray(upperb,&x_ptr);
274: x_ptr[0] = 1.1;
275: VecRestoreArray(upperb,&x_ptr);
276: TaoSetVariableBounds(tao,lowerb,upperb);
278: /* Check for any TAO command line options */
279: TaoSetFromOptions(tao);
280: TaoGetKSP(tao,&ksp);
281: if (ksp) {
282: KSPGetPC(ksp,&pc);
283: PCSetType(pc,PCNONE);
284: }
286: TaoSetTolerances(tao,1e-15,1e-15,1e-15,1e-15,1e-15);
287: /* SOLVE THE APPLICATION */
288: TaoSolve(tao);
290: /* Get information on termination */
291: TaoGetConvergedReason(tao,&reason);
292: if (reason <= 0){
293: ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
294: }
296: VecView(p,PETSC_VIEWER_STDOUT_WORLD);
297: VecDestroy(&p);
298: VecDestroy(&lowerb);
299: VecDestroy(&upperb);
300: TaoDestroy(&tao);
301: PetscFinalize();
302: return 0;
303: }
305: /* ------------------------------------------------------------------ */
308: /*
309: FormFunctionGradient - Evaluates the function and corresponding gradient.
311: Input Parameters:
312: tao - the Tao context
313: X - the input vector
314: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
316: Output Parameters:
317: f - the newly evaluated function
318: G - the newly evaluated gradient
319: */
320: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
321: {
322: AppCtx *ctx = (AppCtx*)ctx0;
323: TS ts;
324: Vec U; /* solution will be stored here */
325: Mat A; /* Jacobian matrix */
326: Mat Jacp; /* Jacobian matrix */
328: PetscInt n = 2;
329: PetscReal ftime;
330: PetscInt steps;
331: PetscScalar *u;
332: PetscScalar *x_ptr,*y_ptr;
333: Vec lambda[1],q,mu[1];
335: VecGetArray(P,&x_ptr);
336: ctx->Pm = x_ptr[0];
337: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
338: Create necessary matrix and vectors
339: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
340: MatCreate(PETSC_COMM_WORLD,&A);
341: MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
342: MatSetType(A,MATDENSE);
343: MatSetFromOptions(A);
344: MatSetUp(A);
346: MatCreateVecs(A,&U,NULL);
348: MatCreate(PETSC_COMM_WORLD,&Jacp);
349: MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
350: MatSetFromOptions(Jacp);
351: MatSetUp(Jacp);
353: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
354: Create timestepping solver context
355: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
356: TSCreate(PETSC_COMM_WORLD,&ts);
357: TSSetProblemType(ts,TS_NONLINEAR);
358: TSSetType(ts,TSCN);
359: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
360: TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,ctx);
361: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
363: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
364: Set initial conditions
365: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
366: VecGetArray(U,&u);
367: u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
368: u[1] = 1.0;
369: VecRestoreArray(U,&u);
370: TSSetSolution(ts,U);
372: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
373: Save trajectory of solution so that TSAdjointSolve() may be used
374: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
375: TSSetSaveTrajectory(ts);
377: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
378: Set solver options
379: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
380: TSSetDuration(ts,PETSC_DEFAULT,10.0);
381: TSSetInitialTimeStep(ts,0.0,.01);
382: TSSetFromOptions(ts);
384: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
385: Solve nonlinear system
386: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
387: TSSolve(ts,U);
389: TSGetSolveTime(ts,&ftime);
390: TSGetTimeStepNumber(ts,&steps);
391: /* VecView(U,PETSC_VIEWER_STDOUT_WORLD); */
393: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
394: Adjoint model starts here
395: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
396: MatCreateVecs(A,&lambda[0],NULL);
397: /* Set initial conditions for the adjoint integration */
398: VecGetArray(lambda[0],&y_ptr);
399: y_ptr[0] = 0.0; y_ptr[1] = 0.0;
400: VecRestoreArray(lambda[0],&y_ptr);
402: MatCreateVecs(Jacp,&mu[0],NULL);
403: VecGetArray(mu[0],&x_ptr);
404: x_ptr[0] = -1.0;
405: VecRestoreArray(mu[0],&x_ptr);
406: TSSetCostGradients(ts,1,lambda,mu);
408: /* Set RHS JacobianP */
409: TSAdjointSetRHSJacobian(ts,Jacp,RHSJacobianP,ctx);
411: TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
412: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
413: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,ctx);
415: TSAdjointSolve(ts);
416: TSGetCostIntegral(ts,&q);
417: /* VecView(q,PETSC_VIEWER_STDOUT_WORLD); */
418: ComputeSensiP(lambda[0],mu[0],ctx);
419: VecCopy(mu[0],G);
421: TSGetCostIntegral(ts,&q);
422: VecGetArray(q,&x_ptr);
423: *f = -ctx->Pm + x_ptr[0];
425: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
426: Free work space. All PETSc objects should be destroyed when they are no longer needed.
427: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
428: MatDestroy(&A);
429: MatDestroy(&Jacp);
430: VecDestroy(&U);
431: VecDestroy(&lambda[0]);
432: VecDestroy(&mu[0]);
433: TSDestroy(&ts);
435: return 0;
436: }