Actual source code: ex3opt.c
petsc-3.8.3 2017-12-09
2: static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\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: This code demonstrates how to solve a ODE-constrained optimization problem with TAO, TSEvent, TSAdjoint and TS.
15: The problem features discontinuities and a cost function in integral form.
16: The gradient is computed with the discrete adjoint of an implicit theta method, see ex3adj.c for details.
17: */
18: #include <petsctao.h>
19: #include <petscts.h>
21: typedef struct {
22: PetscScalar H,D,omega_b,omega_s,Pmax,Pmax_ini,Pm,E,V,X,u_s,c;
23: PetscInt beta;
24: PetscReal tf,tcl;
25: } AppCtx;
27: PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
29: /* Event check */
30: PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx)
31: {
32: AppCtx *user=(AppCtx*)ctx;
35: /* Event for fault-on time */
36: fvalue[0] = t - user->tf;
37: /* Event for fault-off time */
38: fvalue[1] = t - user->tcl;
40: return(0);
41: }
43: PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx)
44: {
45: AppCtx *user=(AppCtx*)ctx;
49: if (event_list[0] == 0) {
50: if (forwardsolve) user->Pmax = 0.0; /* Apply disturbance - this is done by setting Pmax = 0 */
51: else user->Pmax = user->Pmax_ini; /* Going backward, reversal of event */
52: } else if(event_list[0] == 1) {
53: if (forwardsolve) user->Pmax = user->Pmax_ini; /* Remove the fault - this is done by setting Pmax = Pmax_ini */
54: else user->Pmax = 0.0; /* Going backward, reversal of event */
55: }
56: return(0);
57: }
59: /*
60: Defines the ODE passed to the ODE solver
61: */
62: static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx)
63: {
64: PetscErrorCode ierr;
65: PetscScalar *f,Pmax;
66: const PetscScalar *u,*udot;
69: /* The next three lines allow us to access the entries of the vectors directly */
70: VecGetArrayRead(U,&u);
71: VecGetArrayRead(Udot,&udot);
72: VecGetArray(F,&f);
73: Pmax = ctx->Pmax;
74: f[0] = udot[0] - ctx->omega_b*(u[1] - ctx->omega_s);
75: f[1] = 2.0*ctx->H/ctx->omega_s*udot[1] + Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - ctx->omega_s)- ctx->Pm;
77: VecRestoreArrayRead(U,&u);
78: VecRestoreArrayRead(Udot,&udot);
79: VecRestoreArray(F,&f);
80: return(0);
81: }
83: /*
84: Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian.
85: */
86: static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx)
87: {
88: PetscErrorCode ierr;
89: PetscInt rowcol[] = {0,1};
90: PetscScalar J[2][2],Pmax;
91: const PetscScalar *u,*udot;
94: VecGetArrayRead(U,&u);
95: VecGetArrayRead(Udot,&udot);
96: Pmax = ctx->Pmax;
97: J[0][0] = a; J[0][1] = -ctx->omega_b;
98: J[1][1] = 2.0*ctx->H/ctx->omega_s*a + ctx->D; J[1][0] = Pmax*PetscCosScalar(u[0]);
100: MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);
101: VecRestoreArrayRead(U,&u);
102: VecRestoreArrayRead(Udot,&udot);
104: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
105: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
106: if (A != B) {
107: MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);
108: MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);
109: }
110: return(0);
111: }
113: static PetscErrorCode RHSJacobianP(TS ts,PetscReal t,Vec X,Mat A,void *ctx0)
114: {
116: PetscInt row[] = {0,1},col[]={0};
117: PetscScalar J[2][1];
120: J[0][0] = 0;
121: J[1][0] = 1.;
122: MatSetValues(A,2,row,1,col,&J[0][0],INSERT_VALUES);
123: MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);
124: MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);
125: return(0);
126: }
128: static PetscErrorCode CostIntegrand(TS ts,PetscReal t,Vec U,Vec R,AppCtx *ctx)
129: {
130: PetscErrorCode ierr;
131: PetscScalar *r;
132: const PetscScalar *u;
135: VecGetArrayRead(U,&u);
136: VecGetArray(R,&r);
137: r[0] = ctx->c*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta);
138: VecRestoreArray(R,&r);
139: VecRestoreArrayRead(U,&u);
140: return(0);
141: }
143: static PetscErrorCode DRDYFunction(TS ts,PetscReal t,Vec U,Vec *drdy,AppCtx *ctx)
144: {
145: PetscErrorCode ierr;
146: PetscScalar *ry;
147: const PetscScalar *u;
150: VecGetArrayRead(U,&u);
151: VecGetArray(drdy[0],&ry);
152: ry[0] = ctx->c*ctx->beta*PetscPowScalarInt(PetscMax(0., u[0]-ctx->u_s),ctx->beta-1);
153: VecRestoreArray(drdy[0],&ry);
154: VecRestoreArrayRead(U,&u);
155: return(0);
156: }
158: static 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: }
170: PetscErrorCode ComputeSensiP(Vec lambda,Vec mu,AppCtx *ctx)
171: {
172: PetscErrorCode ierr;
173: PetscScalar *y,sensip;
174: const PetscScalar *x;
177: VecGetArrayRead(lambda,&x);
178: VecGetArray(mu,&y);
179: sensip = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax*x[0]+y[0];
180: /* PetscPrintf(PETSC_COMM_WORLD,"\n sensitivity wrt parameter pm: %g \n",(double)sensip); */
181: y[0] = sensip;
182: VecRestoreArray(mu,&y);
183: VecRestoreArrayRead(lambda,&x);
184: return(0);
185: }
187: PetscErrorCode monitor(Tao tao,AppCtx *ctx)
188: {
189: FILE *fp;
190: PetscInt iterate;
191: PetscReal f,gnorm,cnorm,xdiff;
192: TaoConvergedReason reason;
193: PetscErrorCode ierr;
196: TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason);
198: fp = fopen("ex3opt_conv.out","a");
199: PetscFPrintf(PETSC_COMM_WORLD,fp,"%D %g\n",iterate,(double)gnorm);
200: fclose(fp);
201: return(0);
202: }
204: int main(int argc,char **argv)
205: {
206: Vec p;
207: PetscScalar *x_ptr;
208: PetscErrorCode ierr;
209: PetscMPIInt size;
210: AppCtx ctx;
211: Tao tao;
212: KSP ksp;
213: PC pc;
214: Vec lowerb,upperb;
215: PetscBool printtofile;
217: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
218: Initialize program
219: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
220: PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
222: MPI_Comm_size(PETSC_COMM_WORLD,&size);
223: if (size != 1) SETERRQ(PETSC_COMM_SELF,1,"This is a uniprocessor example only!");
225: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226: Set runtime options
227: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
228: PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");
229: {
230: ctx.beta = 2;
231: ctx.c = 10000.0;
232: ctx.u_s = 1.0;
233: ctx.omega_s = 1.0;
234: ctx.omega_b = 120.0*PETSC_PI;
235: ctx.H = 5.0;
236: PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);
237: ctx.D = 5.0;
238: PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);
239: ctx.E = 1.1378;
240: ctx.V = 1.0;
241: ctx.X = 0.545;
242: ctx.Pmax = ctx.E*ctx.V/ctx.X;;
243: ctx.Pmax_ini = ctx.Pmax;
244: PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);
245: ctx.Pm = 1.06;
246: PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);
247: ctx.tf = 0.1;
248: ctx.tcl = 0.2;
249: PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);
250: PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);
251: printtofile = PETSC_FALSE;
252: PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL);
253: }
254: PetscOptionsEnd();
256: /* Create TAO solver and set desired solution method */
257: TaoCreate(PETSC_COMM_WORLD,&tao);
258: TaoSetType(tao,TAOBLMVM);
259: if(printtofile) {
260: TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL);
261: }
262: /*
263: Optimization starts
264: */
265: /* Set initial solution guess */
266: VecCreateSeq(PETSC_COMM_WORLD,1,&p);
267: VecGetArray(p,&x_ptr);
268: x_ptr[0] = ctx.Pm;
269: VecRestoreArray(p,&x_ptr);
271: TaoSetInitialVector(tao,p);
272: /* Set routine for function and gradient evaluation */
273: TaoSetObjectiveAndGradientRoutine(tao,FormFunctionGradient,(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: x_ptr[0] = 1.1;
283: VecRestoreArray(upperb,&x_ptr);
284: TaoSetVariableBounds(tao,lowerb,upperb);
286: /* Check for any TAO command line options */
287: TaoSetFromOptions(tao);
288: TaoGetKSP(tao,&ksp);
289: if (ksp) {
290: KSPGetPC(ksp,&pc);
291: PCSetType(pc,PCNONE);
292: }
294: /* SOLVE THE APPLICATION */
295: TaoSolve(tao);
297: VecView(p,PETSC_VIEWER_STDOUT_WORLD);
298: VecDestroy(&p);
299: VecDestroy(&lowerb);
300: VecDestroy(&upperb);
301: TaoDestroy(&tao);
302: PetscFinalize();
303: return ierr;
304: }
306: /* ------------------------------------------------------------------ */
307: /*
308: FormFunctionGradient - Evaluates the function and corresponding gradient.
310: Input Parameters:
311: tao - the Tao context
312: X - the input vector
313: ptr - optional user-defined context, as set by TaoSetObjectiveAndGradientRoutine()
315: Output Parameters:
316: f - the newly evaluated function
317: G - the newly evaluated gradient
318: */
319: PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
320: {
321: AppCtx *ctx = (AppCtx*)ctx0;
322: TS ts;
323: Vec U; /* solution will be stored here */
324: Mat A; /* Jacobian matrix */
325: Mat Jacp; /* Jacobian matrix */
327: PetscInt n = 2;
328: PetscReal ftime;
329: PetscInt steps;
330: PetscScalar *u;
331: PetscScalar *x_ptr,*y_ptr;
332: Vec lambda[1],q,mu[1];
333: PetscInt direction[2];
334: PetscBool terminate[2];
336: VecGetArray(P,&x_ptr);
337: ctx->Pm = x_ptr[0];
338: VecRestoreArray(P,&x_ptr);
340: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
341: Create necessary matrix and vectors
342: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
343: MatCreate(PETSC_COMM_WORLD,&A);
344: MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);
345: MatSetType(A,MATDENSE);
346: MatSetFromOptions(A);
347: MatSetUp(A);
349: MatCreateVecs(A,&U,NULL);
351: MatCreate(PETSC_COMM_WORLD,&Jacp);
352: MatSetSizes(Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1);
353: MatSetFromOptions(Jacp);
354: MatSetUp(Jacp);
356: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
357: Create timestepping solver context
358: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
359: TSCreate(PETSC_COMM_WORLD,&ts);
360: TSSetProblemType(ts,TS_NONLINEAR);
361: TSSetType(ts,TSCN);
362: TSSetIFunction(ts,NULL,(TSIFunction) IFunction,ctx);
363: TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,ctx);
365: TSSetCostIntegrand(ts,1,NULL,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
366: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
367: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_TRUE,ctx);
369: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
370: Set initial conditions
371: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
372: VecGetArray(U,&u);
373: u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
374: u[1] = 1.0;
375: VecRestoreArray(U,&u);
376: TSSetSolution(ts,U);
378: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
379: Save trajectory of solution so that TSAdjointSolve() may be used
380: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
381: TSSetSaveTrajectory(ts);
383: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
384: Set solver options
385: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
386: TSSetMaxTime(ts,1.0);
387: TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);
388: TSSetTimeStep(ts,0.03125);
389: TSSetFromOptions(ts);
391: direction[0] = direction[1] = 1;
392: terminate[0] = terminate[1] = PETSC_FALSE;
394: TSSetEventHandler(ts,2,direction,terminate,EventFunction,PostEventFunction,(void*)ctx);
396: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
397: Solve nonlinear system
398: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
399: TSSolve(ts,U);
401: TSGetSolveTime(ts,&ftime);
402: TSGetStepNumber(ts,&steps);
403: /* VecView(U,PETSC_VIEWER_STDOUT_WORLD); */
405: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
406: Adjoint model starts here
407: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
408: MatCreateVecs(A,&lambda[0],NULL);
409: /* Set initial conditions for the adjoint integration */
410: VecGetArray(lambda[0],&y_ptr);
411: y_ptr[0] = 0.0; y_ptr[1] = 0.0;
412: VecRestoreArray(lambda[0],&y_ptr);
414: MatCreateVecs(Jacp,&mu[0],NULL);
415: VecGetArray(mu[0],&x_ptr);
416: x_ptr[0] = -1.0;
417: VecRestoreArray(mu[0],&x_ptr);
418: TSSetCostGradients(ts,1,lambda,mu);
420: /* Set RHS JacobianP */
421: TSAdjointSetRHSJacobian(ts,Jacp,RHSJacobianP,ctx);
423: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
424: One can set up the integral to be evaluated during the forward run
425: instead by calling this function before TSSolve and specifying
426: PETSC_TRUE for the second last argument
427: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
428: TSSetCostIntegrand(ts,1,NULL,(PetscErrorCode (*)(TS,PetscReal,Vec,Vec,void*))CostIntegrand,
429: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDYFunction,
430: (PetscErrorCode (*)(TS,PetscReal,Vec,Vec*,void*))DRDPFunction,PETSC_FALSE,ctx);
432: TSAdjointSolve(ts);
433: TSGetCostIntegral(ts,&q);
434: ComputeSensiP(lambda[0],mu[0],ctx);
435: VecCopy(mu[0],G);
437: TSGetCostIntegral(ts,&q);
438: VecGetArray(q,&x_ptr);
439: *f = -ctx->Pm + x_ptr[0];
440: VecRestoreArray(q,&x_ptr);
442: /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
443: Free work space. All PETSc objects should be destroyed when they are no longer needed.
444: - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
445: MatDestroy(&A);
446: MatDestroy(&Jacp);
447: VecDestroy(&U);
448: VecDestroy(&lambda[0]);
449: VecDestroy(&mu[0]);
450: TSDestroy(&ts);
452: return 0;
453: }