Actual source code: ex3opt_fd.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 FormFunction(Tao,Vec,PetscReal*,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 ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
163: {
164: PetscErrorCode ierr;
165: PetscScalar *y,sensip;
166: const PetscScalar *x;
169: VecGetArrayRead(lambda,&x);
170: VecGetArray(mu,&y);
171: sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
172: /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %g \n",(double)sensip); */
173: y[0] = sensip;
174: VecRestoreArray(mu,&y);
175: VecRestoreArrayRead(lambda,&x);
176: return(0);
177: }
181: int main(int argc,char **argv)
182: {
183: Vec p;
184: PetscScalar *x_ptr;
185: PetscErrorCode ierr;
186: PetscMPIInt size;
187: AppCtx ctx;
188: Vec lowerb,upperb;
189: Tao tao;
190: TaoConvergedReason reason;
191: KSP ksp;
192: PC pc;
194: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
195: Initialize program
196: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
197: PetscInitialize(&argc,&argv,NULL,help);
199: MPI_Comm_size(PETSC_COMM_WORLD,&size);
200: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
202: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
203: Set runtime options
204: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
205: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
206: {
207: ctx.beta = 2;
208: ctx.c = 10000.0;
209: ctx.u_s = 1.0;
210: ctx.omega_s = 1.0;
211: ctx.omega_b = 120.0*PETSC_PI;
212: ctx.H = 5.0;
213: PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
214: ctx.D = 5.0;
215: PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
216: ctx.E = 1.1378;
217: ctx.V = 1.0;
218: ctx.X = 0.545;
219: ctx.Pmax = ctx.E*ctx.V/ctx.X;;
220: PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
221: ctx.Pm = 1.061605;
222: PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
223: ctx.tf = 0.1;
224: ctx.tcl = 0.2;
225: PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
226: PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
228: }
229: PetscOptionsEnd();
231: /* Create TAO solver and set desired solution method */
232: TaoCreate(PETSC_COMM_WORLD,&tao);
233: TaoSetType(tao,TAOBLMVM);
235: /*
236: Optimization starts
237: */
238: /* Set initial solution guess */
239: VecCreateSeq(PETSC_COMM_WORLD,1,&p);
240: VecGetArray(p,&x_ptr);
241: x_ptr[0] = ctx.Pm;
242: VecRestoreArray(p,&x_ptr);
244: TaoSetInitialVector(tao,p);
245: /* Set routine for function and gradient evaluation */
246: TaoSetObjectiveRoutine(tao,FormFunction,(void *)&ctx);
247: TaoSetGradientRoutine(tao,TaoDefaultComputeGradient,(void *)&ctx);
249: /* Set bounds for the optimization */
250: VecDuplicate(p,&lowerb);
251: VecDuplicate(p,&upperb);
252: VecGetArray(lowerb,&x_ptr);
253: x_ptr[0] = 0.;
254: VecRestoreArray(lowerb,&x_ptr);
255: VecGetArray(upperb,&x_ptr);
256: x_ptr[0] = 1.1;;
257: VecRestoreArray(upperb,&x_ptr);
258: TaoSetVariableBounds(tao,lowerb,upperb);
260: /* Check for any TAO command line options */
261: TaoSetFromOptions(tao);
262: TaoGetKSP(tao,&ksp);
263: if (ksp) {
264: KSPGetPC(ksp,&pc);
265: PCSetType(pc,PCNONE);
266: }
268: TaoSetTolerances(tao,1e-15,1e-15,1e-15,1e-15,1e-15);
269: /* SOLVE THE APPLICATION */
270: TaoSolve(tao);
272: /* Get information on termination */
273: TaoGetConvergedReason(tao,&reason);
274: if (reason <= 0){
275: ierr=PetscPrintf(MPI_COMM_WORLD, "Try another method! \n");
276: }
278: VecView(p,PETSC_VIEWER_STDOUT_WORLD);
279:
280: /* Free TAO data structures */
281: TaoDestroy(&tao);
282: VecDestroy(&p);
283: VecDestroy(&lowerb);
284: VecDestroy(&upperb);
285: PetscFinalize();
286: return 0;
287: }
289: /* ------------------------------------------------------------------ */
292: /*
293: FormFunction - Evaluates the function and corresponding gradient.
295: Input Parameters:
296: tao - the Tao context
297: X - the input vector
298: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
300: Output Parameters:
301: f - the newly evaluated function
302: */
303: PetscErrorCode FormFunction(Tao tao,Vec P,PetscReal *f,void *ctx0)
304: {
305: AppCtx *ctx = (AppCtx*)ctx0;
306: TS ts;
308: Vec U; /* solution will be stored here */
309: Mat A; /* Jacobian matrix */
310: Mat Jacp; /* Jacobian matrix */
312: PetscInt n = 2;
313: PetscReal ftime;
314: PetscInt steps;
315: PetscScalar *u;
316: PetscScalar *x_ptr,*y_ptr;
317: Vec lambda[1],q,mu[1];
319: VecGetArray(P,&x_ptr);
320: ctx->Pm = x_ptr[0];
321: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
322: Create necessary matrix and vectors
323: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
324: MatCreate(PETSC_COMM_WORLD,&A);
325: MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
326: MatSetType(A,MATDENSE);
327: MatSetFromOptions(A);
328: MatSetUp(A);
330: MatCreateVecs(A,&U,NULL);
332: MatCreate(PETSC_COMM_WORLD,&Jacp);
333: MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
334: MatSetFromOptions(Jacp);
335: MatSetUp(Jacp);
337: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
338: Create timestepping solver context
339: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
340: TSCreate(PETSC_COMM_WORLD,&ts);
341: TSSetProblemType(ts,TS_NONLINEAR);
342: TSSetType(ts,TSCN);
343: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
344: TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,ctx);
345: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
347: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
348: Set initial conditions
349: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
350: VecGetArray(U,&u);
351: u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
352: u[1] = 1.0;
353: VecRestoreArray(U,&u);
354: TSSetSolution(ts,U);
356: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
357: Save trajectory of solution so that TSAdjointSolve() may be used
358: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
359: TSSetSaveTrajectory(ts);
361: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
362: Set solver options
363: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
364: TSSetDuration(ts,PETSC_DEFAULT,10.0);
365: TSSetInitialTimeStep(ts,0.0,.01);
366: TSSetFromOptions(ts);
368: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
369: Solve nonlinear system
370: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
371: TSSolve(ts,U);
373: TSGetSolveTime(ts,&ftime);
374: TSGetTimeStepNumber(ts,&steps);
376: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
377: Adjoint model starts here
378: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
379: MatCreateVecs(A,&lambda[0],NULL);
380: /* Set initial conditions for the adjoint integration */
381: VecGetArray(lambda[0],&y_ptr);
382: y_ptr[0] = 0.0; y_ptr[1] = 0.0;
383: VecRestoreArray(lambda[0],&y_ptr);
385: MatCreateVecs(Jacp,&mu[0],NULL);
386: VecGetArray(mu[0],&x_ptr);
387: x_ptr[0] = -1.0;
388: VecRestoreArray(mu[0],&x_ptr);
389: TSSetCostGradients(ts,1,lambda,mu);
391: TSAdjointSetRHSJacobian(ts,Jacp,RHSJacobianP,ctx);
393: TSSetCostIntegrand(ts,1,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
394: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
395: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,ctx);
397: TSAdjointSolve(ts);
398: TSGetCostIntegral(ts,&q);
399: ComputeSensiP(lambda[0],mu[0],ctx);
401: VecGetArray(q,&x_ptr);
402: *f = -ctx->Pm + x_ptr[0];
404: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
405: Free work space. All PETSc objects should be destroyed when they are no longer needed.
406: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
407: MatDestroy(&A);
408: MatDestroy(&Jacp);
409: VecDestroy(&U);
410: VecDestroy(&lambda[0]);
411: VecDestroy(&mu[0]);
412: TSDestroy(&ts);
414: return 0;
415: }