xref: /petsc/src/ts/tutorials/power_grid/ex3opt_fd.c (revision fbf9dbe564678ed6eff1806adbc4c4f01b9743f4)
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   Solve the same optimization problem as in ex3opt.c.
15   Use finite difference to approximate the gradients.
16 */
17 #include <petsctao.h>
18 #include <petscts.h>
19 #include "ex3.h"
20 
21 PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
22 
23 PetscErrorCode monitor(Tao tao, AppCtx *ctx)
24 {
25   FILE              *fp;
26   PetscInt           iterate;
27   PetscReal          f, gnorm, cnorm, xdiff;
28   Vec                X, G;
29   const PetscScalar *x, *g;
30   TaoConvergedReason reason;
31 
32   PetscFunctionBeginUser;
33   PetscCall(TaoGetSolutionStatus(tao, &iterate, &f, &gnorm, &cnorm, &xdiff, &reason));
34   PetscCall(TaoGetSolution(tao, &X));
35   PetscCall(TaoGetGradient(tao, &G, NULL, NULL));
36   PetscCall(VecGetArrayRead(X, &x));
37   PetscCall(VecGetArrayRead(G, &g));
38   fp = fopen("ex3opt_fd_conv.out", "a");
39   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fp, "%" PetscInt_FMT " %g %.12lf %.12lf\n", iterate, (double)gnorm, (double)PetscRealPart(x[0]), (double)PetscRealPart(g[0])));
40   PetscCall(VecRestoreArrayRead(X, &x));
41   PetscCall(VecRestoreArrayRead(G, &g));
42   fclose(fp);
43   PetscFunctionReturn(PETSC_SUCCESS);
44 }
45 
46 int main(int argc, char **argv)
47 {
48   Vec          p;
49   PetscScalar *x_ptr;
50   PetscMPIInt  size;
51   AppCtx       ctx;
52   Vec          lowerb, upperb;
53   Tao          tao;
54   KSP          ksp;
55   PC           pc;
56   PetscBool    printtofile;
57   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58      Initialize program
59      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60   PetscFunctionBeginUser;
61   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
62   PetscFunctionBeginUser;
63   PetscCallMPI(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   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
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     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
78     ctx.D = 5.0;
79     PetscCall(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     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
86     ctx.Pm = 1.06;
87     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
88     ctx.tf  = 0.1;
89     ctx.tcl = 0.2;
90     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
91     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
92     printtofile = PETSC_FALSE;
93     PetscCall(PetscOptionsBool("-printtofile", "Print convergence results to file", "", printtofile, &printtofile, NULL));
94   }
95   PetscOptionsEnd();
96 
97   /* Create TAO solver and set desired solution method */
98   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
99   PetscCall(TaoSetType(tao, TAOBLMVM));
100   if (printtofile) PetscCall(TaoSetMonitor(tao, (PetscErrorCode(*)(Tao, void *))monitor, (void *)&ctx, NULL));
101   PetscCall(TaoSetMaximumIterations(tao, 30));
102   /*
103      Optimization starts
104   */
105   /* Set initial solution guess */
106   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p));
107   PetscCall(VecGetArray(p, &x_ptr));
108   x_ptr[0] = ctx.Pm;
109   PetscCall(VecRestoreArray(p, &x_ptr));
110 
111   PetscCall(TaoSetSolution(tao, p));
112   /* Set routine for function and gradient evaluation */
113   PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx));
114   PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, (void *)&ctx));
115 
116   /* Set bounds for the optimization */
117   PetscCall(VecDuplicate(p, &lowerb));
118   PetscCall(VecDuplicate(p, &upperb));
119   PetscCall(VecGetArray(lowerb, &x_ptr));
120   x_ptr[0] = 0.;
121   PetscCall(VecRestoreArray(lowerb, &x_ptr));
122   PetscCall(VecGetArray(upperb, &x_ptr));
123   x_ptr[0] = 1.1;
124   PetscCall(VecRestoreArray(upperb, &x_ptr));
125   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));
126 
127   /* Check for any TAO command line options */
128   PetscCall(TaoSetFromOptions(tao));
129   PetscCall(TaoGetKSP(tao, &ksp));
130   if (ksp) {
131     PetscCall(KSPGetPC(ksp, &pc));
132     PetscCall(PCSetType(pc, PCNONE));
133   }
134 
135   /* SOLVE THE APPLICATION */
136   PetscCall(TaoSolve(tao));
137 
138   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
139 
140   /* Free TAO data structures */
141   PetscCall(TaoDestroy(&tao));
142   PetscCall(VecDestroy(&p));
143   PetscCall(VecDestroy(&lowerb));
144   PetscCall(VecDestroy(&upperb));
145   PetscCall(PetscFinalize());
146   return 0;
147 }
148 
149 /* ------------------------------------------------------------------ */
150 /*
151    FormFunction - Evaluates the function and corresponding gradient.
152 
153    Input Parameters:
154    tao - the Tao context
155    X   - the input vector
156    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
157 
158    Output Parameters:
159    f   - the newly evaluated function
160 */
161 PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0)
162 {
163   AppCtx            *ctx = (AppCtx *)ctx0;
164   TS                 ts, quadts;
165   Vec                U; /* solution will be stored here */
166   Mat                A; /* Jacobian matrix */
167   PetscInt           n = 2;
168   PetscReal          ftime;
169   PetscInt           steps;
170   PetscScalar       *u;
171   const PetscScalar *x_ptr, *qx_ptr;
172   Vec                q;
173   PetscInt           direction[2];
174   PetscBool          terminate[2];
175 
176   PetscFunctionBeginUser;
177   PetscCall(VecGetArrayRead(P, &x_ptr));
178   ctx->Pm = x_ptr[0];
179   PetscCall(VecRestoreArrayRead(P, &x_ptr));
180   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181     Create necessary matrix and vectors
182     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
183   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
184   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
185   PetscCall(MatSetType(A, MATDENSE));
186   PetscCall(MatSetFromOptions(A));
187   PetscCall(MatSetUp(A));
188 
189   PetscCall(MatCreateVecs(A, &U, NULL));
190 
191   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192      Create timestepping solver context
193      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
194   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
195   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
196   PetscCall(TSSetType(ts, TSCN));
197   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, ctx));
198   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, ctx));
199 
200   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201      Set initial conditions
202    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
203   PetscCall(VecGetArray(U, &u));
204   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
205   u[1] = 1.0;
206   PetscCall(VecRestoreArray(U, &u));
207   PetscCall(TSSetSolution(ts, U));
208 
209   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210      Set solver options
211    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
212   PetscCall(TSSetMaxTime(ts, 1.0));
213   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
214   PetscCall(TSSetTimeStep(ts, 0.03125));
215   PetscCall(TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts));
216   PetscCall(TSGetSolution(quadts, &q));
217   PetscCall(VecSet(q, 0.0));
218   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, ctx));
219   PetscCall(TSSetFromOptions(ts));
220 
221   direction[0] = direction[1] = 1;
222   terminate[0] = terminate[1] = PETSC_FALSE;
223 
224   PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)ctx));
225 
226   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227      Solve nonlinear system
228      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
229   PetscCall(TSSolve(ts, U));
230 
231   PetscCall(TSGetSolveTime(ts, &ftime));
232   PetscCall(TSGetStepNumber(ts, &steps));
233   PetscCall(VecGetArrayRead(q, &qx_ptr));
234   *f = -ctx->Pm + qx_ptr[0];
235   PetscCall(VecRestoreArrayRead(q, &qx_ptr));
236 
237   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
239    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
240   PetscCall(MatDestroy(&A));
241   PetscCall(VecDestroy(&U));
242   PetscCall(TSDestroy(&ts));
243   PetscFunctionReturn(PETSC_SUCCESS);
244 }
245 
246 /*TEST
247 
248    build:
249       requires: !complex !single
250 
251    test:
252       args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3
253 
254 TEST*/
255