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