xref: /petsc/src/ts/tutorials/power_grid/ex3opt.c (revision 40badf4fbc550ac1f60bd080eaff6de6d55b946d)
1 
2 static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\n";
3 
4 /*F
5 
6 \begin{eqnarray}
7                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9 \end{eqnarray}
10 
11 F*/
12 
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 
19 #include <petsctao.h>
20 #include <petscts.h>
21 #include "ex3.h"
22 
23 PetscErrorCode FormFunctionGradient(Tao,Vec,PetscReal*,Vec,void*);
24 
25 PetscErrorCode monitor(Tao tao,AppCtx *ctx)
26 {
27   FILE               *fp;
28   PetscInt           iterate;
29   PetscReal          f,gnorm,cnorm,xdiff;
30   TaoConvergedReason reason;
31 
32   PetscFunctionBeginUser;
33   CHKERRQ(TaoGetSolutionStatus(tao,&iterate,&f,&gnorm,&cnorm,&xdiff,&reason));
34 
35   fp = fopen("ex3opt_conv.out","a");
36   CHKERRQ(PetscFPrintf(PETSC_COMM_WORLD,fp,"%D %g\n",iterate,(double)gnorm));
37   fclose(fp);
38   PetscFunctionReturn(0);
39 }
40 
41 int main(int argc,char **argv)
42 {
43   Vec                p;
44   PetscScalar        *x_ptr;
45   PetscErrorCode     ierr;
46   PetscMPIInt        size;
47   AppCtx             ctx;
48   Tao                tao;
49   KSP                ksp;
50   PC                 pc;
51   Vec                lambda[1],mu[1],lowerb,upperb;
52   PetscBool          printtofile;
53   PetscInt           direction[2];
54   PetscBool          terminate[2];
55   Mat                qgrad;         /* Forward sesivitiy */
56   Mat                sp;            /* Forward sensitivity matrix */
57 
58   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
59      Initialize program
60      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
61   ierr = PetscInitialize(&argc,&argv,NULL,help);if (ierr) return ierr;
62   PetscFunctionBeginUser;
63   CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size));
64   PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"This is a uniprocessor example only!");
65 
66   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67     Set runtime options
68     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69   ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr);
70   {
71     ctx.beta    = 2;
72     ctx.c       = 10000.0;
73     ctx.u_s     = 1.0;
74     ctx.omega_s = 1.0;
75     ctx.omega_b = 120.0*PETSC_PI;
76     ctx.H       = 5.0;
77     CHKERRQ(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL));
78     ctx.D       = 5.0;
79     CHKERRQ(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL));
80     ctx.E       = 1.1378;
81     ctx.V       = 1.0;
82     ctx.X       = 0.545;
83     ctx.Pmax    = ctx.E*ctx.V/ctx.X;
84     ctx.Pmax_ini = ctx.Pmax;
85     CHKERRQ(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL));
86     ctx.Pm      = 1.06;
87     CHKERRQ(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL));
88     ctx.tf      = 0.1;
89     ctx.tcl     = 0.2;
90     CHKERRQ(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL));
91     CHKERRQ(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL));
92     printtofile = PETSC_FALSE;
93     CHKERRQ(PetscOptionsBool("-printtofile","Print convergence results to file","",printtofile,&printtofile,NULL));
94     ctx.sa      = SA_ADJ;
95     CHKERRQ(PetscOptionsEnum("-sa_method","Sensitivity analysis method (adj or tlm)","",SAMethods,(PetscEnum)ctx.sa,(PetscEnum*)&ctx.sa,NULL));
96   }
97   ierr = PetscOptionsEnd();CHKERRQ(ierr);
98 
99   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
100     Create necessary matrix and vectors
101     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
102   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&ctx.Jac));
103   CHKERRQ(MatSetSizes(ctx.Jac,2,2,PETSC_DETERMINE,PETSC_DETERMINE));
104   CHKERRQ(MatSetType(ctx.Jac,MATDENSE));
105   CHKERRQ(MatSetFromOptions(ctx.Jac));
106   CHKERRQ(MatSetUp(ctx.Jac));
107   CHKERRQ(MatCreate(PETSC_COMM_WORLD,&ctx.Jacp));
108   CHKERRQ(MatSetSizes(ctx.Jacp,PETSC_DECIDE,PETSC_DECIDE,2,1));
109   CHKERRQ(MatSetFromOptions(ctx.Jacp));
110   CHKERRQ(MatSetUp(ctx.Jacp));
111   CHKERRQ(MatCreateVecs(ctx.Jac,&ctx.U,NULL));
112   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&ctx.DRDP));
113   CHKERRQ(MatSetUp(ctx.DRDP));
114   CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&ctx.DRDU));
115   CHKERRQ(MatSetUp(ctx.DRDU));
116 
117   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
118      Create timestepping solver context
119      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
120   CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ctx.ts));
121   CHKERRQ(TSSetProblemType(ctx.ts,TS_NONLINEAR));
122   CHKERRQ(TSSetType(ctx.ts,TSCN));
123   CHKERRQ(TSSetRHSFunction(ctx.ts,NULL,(TSRHSFunction)RHSFunction,&ctx));
124   CHKERRQ(TSSetRHSJacobian(ctx.ts,ctx.Jac,ctx.Jac,(TSRHSJacobian)RHSJacobian,&ctx));
125   CHKERRQ(TSSetRHSJacobianP(ctx.ts,ctx.Jacp,RHSJacobianP,&ctx));
126 
127   if (ctx.sa == SA_ADJ) {
128     CHKERRQ(MatCreateVecs(ctx.Jac,&lambda[0],NULL));
129     CHKERRQ(MatCreateVecs(ctx.Jacp,&mu[0],NULL));
130     CHKERRQ(TSSetSaveTrajectory(ctx.ts));
131     CHKERRQ(TSSetCostGradients(ctx.ts,1,lambda,mu));
132     CHKERRQ(TSCreateQuadratureTS(ctx.ts,PETSC_FALSE,&ctx.quadts));
133     CHKERRQ(TSSetRHSFunction(ctx.quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx));
134     CHKERRQ(TSSetRHSJacobian(ctx.quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx));
135     CHKERRQ(TSSetRHSJacobianP(ctx.quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx));
136   }
137   if (ctx.sa == SA_TLM) {
138     CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,1,1,NULL,&qgrad));
139     CHKERRQ(MatCreateDense(PETSC_COMM_WORLD,PETSC_DECIDE,PETSC_DECIDE,2,1,NULL,&sp));
140     CHKERRQ(TSForwardSetSensitivities(ctx.ts,1,sp));
141     CHKERRQ(TSCreateQuadratureTS(ctx.ts,PETSC_TRUE,&ctx.quadts));
142     CHKERRQ(TSForwardSetSensitivities(ctx.quadts,1,qgrad));
143     CHKERRQ(TSSetRHSFunction(ctx.quadts,NULL,(TSRHSFunction)CostIntegrand,&ctx));
144     CHKERRQ(TSSetRHSJacobian(ctx.quadts,ctx.DRDU,ctx.DRDU,(TSRHSJacobian)DRDUJacobianTranspose,&ctx));
145     CHKERRQ(TSSetRHSJacobianP(ctx.quadts,ctx.DRDP,DRDPJacobianTranspose,&ctx));
146   }
147 
148   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
149      Set solver options
150    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
151   CHKERRQ(TSSetMaxTime(ctx.ts,1.0));
152   CHKERRQ(TSSetExactFinalTime(ctx.ts,TS_EXACTFINALTIME_MATCHSTEP));
153   CHKERRQ(TSSetTimeStep(ctx.ts,0.03125));
154   CHKERRQ(TSSetFromOptions(ctx.ts));
155 
156   direction[0] = direction[1] = 1;
157   terminate[0] = terminate[1] = PETSC_FALSE;
158   CHKERRQ(TSSetEventHandler(ctx.ts,2,direction,terminate,EventFunction,PostEventFunction,&ctx));
159 
160   /* Create TAO solver and set desired solution method */
161   CHKERRQ(TaoCreate(PETSC_COMM_WORLD,&tao));
162   CHKERRQ(TaoSetType(tao,TAOBLMVM));
163   if (printtofile) {
164     CHKERRQ(TaoSetMonitor(tao,(PetscErrorCode (*)(Tao, void*))monitor,(void *)&ctx,PETSC_NULL));
165   }
166   /*
167      Optimization starts
168   */
169   /* Set initial solution guess */
170   CHKERRQ(VecCreateSeq(PETSC_COMM_WORLD,1,&p));
171   CHKERRQ(VecGetArray(p,&x_ptr));
172   x_ptr[0] = ctx.Pm;
173   CHKERRQ(VecRestoreArray(p,&x_ptr));
174 
175   CHKERRQ(TaoSetSolution(tao,p));
176   /* Set routine for function and gradient evaluation */
177   CHKERRQ(TaoSetObjectiveAndGradient(tao,NULL,FormFunctionGradient,(void *)&ctx));
178 
179   /* Set bounds for the optimization */
180   CHKERRQ(VecDuplicate(p,&lowerb));
181   CHKERRQ(VecDuplicate(p,&upperb));
182   CHKERRQ(VecGetArray(lowerb,&x_ptr));
183   x_ptr[0] = 0.;
184   CHKERRQ(VecRestoreArray(lowerb,&x_ptr));
185   CHKERRQ(VecGetArray(upperb,&x_ptr));
186   x_ptr[0] = 1.1;
187   CHKERRQ(VecRestoreArray(upperb,&x_ptr));
188   CHKERRQ(TaoSetVariableBounds(tao,lowerb,upperb));
189 
190   /* Check for any TAO command line options */
191   CHKERRQ(TaoSetFromOptions(tao));
192   CHKERRQ(TaoGetKSP(tao,&ksp));
193   if (ksp) {
194     CHKERRQ(KSPGetPC(ksp,&pc));
195     CHKERRQ(PCSetType(pc,PCNONE));
196   }
197 
198   /* SOLVE THE APPLICATION */
199   CHKERRQ(TaoSolve(tao));
200 
201   CHKERRQ(VecView(p,PETSC_VIEWER_STDOUT_WORLD));
202 
203   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
204      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
205    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
206   CHKERRQ(MatDestroy(&ctx.Jac));
207   CHKERRQ(MatDestroy(&ctx.Jacp));
208   CHKERRQ(MatDestroy(&ctx.DRDU));
209   CHKERRQ(MatDestroy(&ctx.DRDP));
210   CHKERRQ(VecDestroy(&ctx.U));
211   if (ctx.sa == SA_ADJ) {
212     CHKERRQ(VecDestroy(&lambda[0]));
213     CHKERRQ(VecDestroy(&mu[0]));
214   }
215   if (ctx.sa == SA_TLM) {
216     CHKERRQ(MatDestroy(&qgrad));
217     CHKERRQ(MatDestroy(&sp));
218   }
219   CHKERRQ(TSDestroy(&ctx.ts));
220   CHKERRQ(VecDestroy(&p));
221   CHKERRQ(VecDestroy(&lowerb));
222   CHKERRQ(VecDestroy(&upperb));
223   CHKERRQ(TaoDestroy(&tao));
224   ierr = PetscFinalize();
225   return ierr;
226 }
227 
228 /* ------------------------------------------------------------------ */
229 /*
230    FormFunctionGradient - Evaluates the function and corresponding gradient.
231 
232    Input Parameters:
233    tao - the Tao context
234    X   - the input vector
235    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
236 
237    Output Parameters:
238    f   - the newly evaluated function
239    G   - the newly evaluated gradient
240 */
241 PetscErrorCode FormFunctionGradient(Tao tao,Vec P,PetscReal *f,Vec G,void *ctx0)
242 {
243   AppCtx         *ctx = (AppCtx*)ctx0;
244   PetscInt       nadj;
245   PetscReal      ftime;
246   PetscInt       steps;
247   PetscScalar    *u;
248   PetscScalar    *x_ptr,*y_ptr;
249   Vec            q;
250   Mat            qgrad;
251 
252   CHKERRQ(VecGetArrayRead(P,(const PetscScalar**)&x_ptr));
253   ctx->Pm = x_ptr[0];
254   CHKERRQ(VecRestoreArrayRead(P,(const PetscScalar**)&x_ptr));
255 
256   /* reinitialize the solution vector */
257   CHKERRQ(VecGetArray(ctx->U,&u));
258   u[0] = PetscAsinScalar(ctx->Pm/ctx->Pmax);
259   u[1] = 1.0;
260   CHKERRQ(VecRestoreArray(ctx->U,&u));
261   CHKERRQ(TSSetSolution(ctx->ts,ctx->U));
262 
263   /* reset time */
264   CHKERRQ(TSSetTime(ctx->ts,0.0));
265 
266   /* reset step counter, this is critical for adjoint solver */
267   CHKERRQ(TSSetStepNumber(ctx->ts,0));
268 
269   /* reset step size, the step size becomes negative after TSAdjointSolve */
270   CHKERRQ(TSSetTimeStep(ctx->ts,0.03125));
271 
272   /* reinitialize the integral value */
273   CHKERRQ(TSGetQuadratureTS(ctx->ts,NULL,&ctx->quadts));
274   CHKERRQ(TSGetSolution(ctx->quadts,&q));
275   CHKERRQ(VecSet(q,0.0));
276 
277   if (ctx->sa == SA_TLM) { /* reset the forward sensitivities */
278     TS             quadts;
279     Mat            sp;
280     PetscScalar    val[2];
281     const PetscInt row[]={0,1},col[]={0};
282 
283     CHKERRQ(TSGetQuadratureTS(ctx->ts,NULL,&quadts));
284     CHKERRQ(TSForwardGetSensitivities(quadts,NULL,&qgrad));
285     CHKERRQ(MatZeroEntries(qgrad));
286     CHKERRQ(TSForwardGetSensitivities(ctx->ts,NULL,&sp));
287     val[0] = 1./PetscSqrtScalar(1.-(ctx->Pm/ctx->Pmax)*(ctx->Pm/ctx->Pmax))/ctx->Pmax;
288     val[1] = 0.0;
289     CHKERRQ(MatSetValues(sp,2,row,1,col,val,INSERT_VALUES));
290     CHKERRQ(MatAssemblyBegin(sp,MAT_FINAL_ASSEMBLY));
291     CHKERRQ(MatAssemblyEnd(sp,MAT_FINAL_ASSEMBLY));
292   }
293 
294   /* solve the ODE */
295   CHKERRQ(TSSolve(ctx->ts,ctx->U));
296   CHKERRQ(TSGetSolveTime(ctx->ts,&ftime));
297   CHKERRQ(TSGetStepNumber(ctx->ts,&steps));
298 
299   if (ctx->sa == SA_ADJ) {
300     Vec *lambda,*mu;
301     /* reset the terminal condition for adjoint */
302     CHKERRQ(TSGetCostGradients(ctx->ts,&nadj,&lambda,&mu));
303     CHKERRQ(VecGetArray(lambda[0],&y_ptr));
304     y_ptr[0] = 0.0; y_ptr[1] = 0.0;
305     CHKERRQ(VecRestoreArray(lambda[0],&y_ptr));
306     CHKERRQ(VecGetArray(mu[0],&x_ptr));
307     x_ptr[0] = -1.0;
308     CHKERRQ(VecRestoreArray(mu[0],&x_ptr));
309 
310     /* solve the adjont */
311     CHKERRQ(TSAdjointSolve(ctx->ts));
312 
313     CHKERRQ(ComputeSensiP(lambda[0],mu[0],ctx));
314     CHKERRQ(VecCopy(mu[0],G));
315   }
316 
317   if (ctx->sa == SA_TLM) {
318     CHKERRQ(VecGetArray(G,&x_ptr));
319     CHKERRQ(MatDenseGetArray(qgrad,&y_ptr));
320     x_ptr[0] = y_ptr[0]-1.;
321     CHKERRQ(MatDenseRestoreArray(qgrad,&y_ptr));
322     CHKERRQ(VecRestoreArray(G,&x_ptr));
323   }
324 
325   CHKERRQ(TSGetSolution(ctx->quadts,&q));
326   CHKERRQ(VecGetArray(q,&x_ptr));
327   *f   = -ctx->Pm + x_ptr[0];
328   CHKERRQ(VecRestoreArray(q,&x_ptr));
329   return 0;
330 }
331 
332 /*TEST
333 
334    build:
335       requires: !complex !single
336 
337    test:
338       args: -viewer_binary_skip_info -ts_type cn -pc_type lu -tao_monitor
339 
340    test:
341       suffix: 2
342       output_file: output/ex3opt_1.out
343       args: -sa_method tlm -ts_type cn -pc_type lu -tao_monitor
344 TEST*/
345