xref: /petsc/src/ts/tutorials/power_grid/ex3opt_fd.c (revision 3ba1676111f5c958fe6c2729b46ca4d523958bb3)
1c4762a1bSJed Brown 
2c4762a1bSJed Brown static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\n";
3c4762a1bSJed Brown 
4c4762a1bSJed Brown /*F
5c4762a1bSJed Brown 
6c4762a1bSJed Brown \begin{eqnarray}
7c4762a1bSJed Brown                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
8c4762a1bSJed Brown                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
9c4762a1bSJed Brown \end{eqnarray}
10c4762a1bSJed Brown 
11c4762a1bSJed Brown F*/
12c4762a1bSJed Brown 
13c4762a1bSJed Brown /*
14c4762a1bSJed Brown   Solve the same optimization problem as in ex3opt.c.
15c4762a1bSJed Brown   Use finite difference to approximate the gradients.
16c4762a1bSJed Brown */
17c4762a1bSJed Brown #include <petsctao.h>
18c4762a1bSJed Brown #include <petscts.h>
19c4762a1bSJed Brown #include "ex3.h"
20c4762a1bSJed Brown 
21c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);
22c4762a1bSJed Brown 
23d71ae5a4SJacob Faibussowitsch PetscErrorCode monitor(Tao tao, AppCtx *ctx)
24d71ae5a4SJacob Faibussowitsch {
25c4762a1bSJed Brown   FILE              *fp;
26c4762a1bSJed Brown   PetscInt           iterate;
27c4762a1bSJed Brown   PetscReal          f, gnorm, cnorm, xdiff;
28c4762a1bSJed Brown   Vec                X, G;
29c4762a1bSJed Brown   const PetscScalar *x, *g;
30c4762a1bSJed Brown   TaoConvergedReason reason;
31c4762a1bSJed Brown 
32c4762a1bSJed Brown   PetscFunctionBeginUser;
339566063dSJacob Faibussowitsch   PetscCall(TaoGetSolutionStatus(tao, &iterate, &f, &gnorm, &cnorm, &xdiff, &reason));
349566063dSJacob Faibussowitsch   PetscCall(TaoGetSolution(tao, &X));
359566063dSJacob Faibussowitsch   PetscCall(TaoGetGradient(tao, &G, NULL, NULL));
369566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(X, &x));
379566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(G, &g));
38c4762a1bSJed Brown   fp = fopen("ex3opt_fd_conv.out", "a");
3963a3b9bcSJacob Faibussowitsch   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fp, "%" PetscInt_FMT " %g %.12lf %.12lf\n", iterate, (double)gnorm, (double)PetscRealPart(x[0]), (double)PetscRealPart(g[0])));
409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(X, &x));
419566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(G, &g));
42c4762a1bSJed Brown   fclose(fp);
43*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
44c4762a1bSJed Brown }
45c4762a1bSJed Brown 
46d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
47d71ae5a4SJacob Faibussowitsch {
48c4762a1bSJed Brown   Vec          p;
49c4762a1bSJed Brown   PetscScalar *x_ptr;
50c4762a1bSJed Brown   PetscMPIInt  size;
51c4762a1bSJed Brown   AppCtx       ctx;
52c4762a1bSJed Brown   Vec          lowerb, upperb;
53c4762a1bSJed Brown   Tao          tao;
54c4762a1bSJed Brown   KSP          ksp;
55c4762a1bSJed Brown   PC           pc;
56c4762a1bSJed Brown   PetscBool    printtofile;
57c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
58c4762a1bSJed Brown      Initialize program
59c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
60327415f7SBarry Smith   PetscFunctionBeginUser;
619566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
62c4762a1bSJed Brown   PetscFunctionBeginUser;
639566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
643c633725SBarry Smith   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");
65c4762a1bSJed Brown 
66c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
67c4762a1bSJed Brown     Set runtime options
68c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
69d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
70c4762a1bSJed Brown   {
71c4762a1bSJed Brown     ctx.beta    = 2;
72c4762a1bSJed Brown     ctx.c       = 10000.0;
73c4762a1bSJed Brown     ctx.u_s     = 1.0;
74c4762a1bSJed Brown     ctx.omega_s = 1.0;
75c4762a1bSJed Brown     ctx.omega_b = 120.0 * PETSC_PI;
76c4762a1bSJed Brown     ctx.H       = 5.0;
779566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
78c4762a1bSJed Brown     ctx.D = 5.0;
799566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
80c4762a1bSJed Brown     ctx.E        = 1.1378;
81c4762a1bSJed Brown     ctx.V        = 1.0;
82c4762a1bSJed Brown     ctx.X        = 0.545;
83c4762a1bSJed Brown     ctx.Pmax     = ctx.E * ctx.V / ctx.X;
84c4762a1bSJed Brown     ctx.Pmax_ini = ctx.Pmax;
859566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
86c4762a1bSJed Brown     ctx.Pm = 1.06;
879566063dSJacob Faibussowitsch     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
88c4762a1bSJed Brown     ctx.tf  = 0.1;
89c4762a1bSJed Brown     ctx.tcl = 0.2;
909566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
919566063dSJacob Faibussowitsch     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
92c4762a1bSJed Brown     printtofile = PETSC_FALSE;
939566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-printtofile", "Print convergence results to file", "", printtofile, &printtofile, NULL));
94c4762a1bSJed Brown   }
95d0609cedSBarry Smith   PetscOptionsEnd();
96c4762a1bSJed Brown 
97c4762a1bSJed Brown   /* Create TAO solver and set desired solution method */
989566063dSJacob Faibussowitsch   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
999566063dSJacob Faibussowitsch   PetscCall(TaoSetType(tao, TAOBLMVM));
100f3fa974cSJacob Faibussowitsch   if (printtofile) PetscCall(TaoSetMonitor(tao, (PetscErrorCode(*)(Tao, void *))monitor, (void *)&ctx, NULL));
1019566063dSJacob Faibussowitsch   PetscCall(TaoSetMaximumIterations(tao, 30));
102c4762a1bSJed Brown   /*
103c4762a1bSJed Brown      Optimization starts
104c4762a1bSJed Brown   */
105c4762a1bSJed Brown   /* Set initial solution guess */
1069566063dSJacob Faibussowitsch   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p));
1079566063dSJacob Faibussowitsch   PetscCall(VecGetArray(p, &x_ptr));
108c4762a1bSJed Brown   x_ptr[0] = ctx.Pm;
1099566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(p, &x_ptr));
110c4762a1bSJed Brown 
1119566063dSJacob Faibussowitsch   PetscCall(TaoSetSolution(tao, p));
112c4762a1bSJed Brown   /* Set routine for function and gradient evaluation */
1139566063dSJacob Faibussowitsch   PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx));
1149566063dSJacob Faibussowitsch   PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, (void *)&ctx));
115c4762a1bSJed Brown 
116c4762a1bSJed Brown   /* Set bounds for the optimization */
1179566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p, &lowerb));
1189566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(p, &upperb));
1199566063dSJacob Faibussowitsch   PetscCall(VecGetArray(lowerb, &x_ptr));
120c4762a1bSJed Brown   x_ptr[0] = 0.;
1219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(lowerb, &x_ptr));
1229566063dSJacob Faibussowitsch   PetscCall(VecGetArray(upperb, &x_ptr));
123c4762a1bSJed Brown   x_ptr[0] = 1.1;
1249566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(upperb, &x_ptr));
1259566063dSJacob Faibussowitsch   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));
126c4762a1bSJed Brown 
127c4762a1bSJed Brown   /* Check for any TAO command line options */
1289566063dSJacob Faibussowitsch   PetscCall(TaoSetFromOptions(tao));
1299566063dSJacob Faibussowitsch   PetscCall(TaoGetKSP(tao, &ksp));
130c4762a1bSJed Brown   if (ksp) {
1319566063dSJacob Faibussowitsch     PetscCall(KSPGetPC(ksp, &pc));
1329566063dSJacob Faibussowitsch     PetscCall(PCSetType(pc, PCNONE));
133c4762a1bSJed Brown   }
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* SOLVE THE APPLICATION */
1369566063dSJacob Faibussowitsch   PetscCall(TaoSolve(tao));
137c4762a1bSJed Brown 
1389566063dSJacob Faibussowitsch   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));
139c4762a1bSJed Brown 
140c4762a1bSJed Brown   /* Free TAO data structures */
1419566063dSJacob Faibussowitsch   PetscCall(TaoDestroy(&tao));
1429566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&p));
1439566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&lowerb));
1449566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&upperb));
1459566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
146b122ec5aSJacob Faibussowitsch   return 0;
147c4762a1bSJed Brown }
148c4762a1bSJed Brown 
149c4762a1bSJed Brown /* ------------------------------------------------------------------ */
150c4762a1bSJed Brown /*
151c4762a1bSJed Brown    FormFunction - Evaluates the function and corresponding gradient.
152c4762a1bSJed Brown 
153c4762a1bSJed Brown    Input Parameters:
154c4762a1bSJed Brown    tao - the Tao context
155c4762a1bSJed Brown    X   - the input vector
156a82e8c82SStefano Zampini    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()
157c4762a1bSJed Brown 
158c4762a1bSJed Brown    Output Parameters:
159c4762a1bSJed Brown    f   - the newly evaluated function
160c4762a1bSJed Brown */
161d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0)
162d71ae5a4SJacob Faibussowitsch {
163c4762a1bSJed Brown   AppCtx            *ctx = (AppCtx *)ctx0;
164c4762a1bSJed Brown   TS                 ts, quadts;
165c4762a1bSJed Brown   Vec                U; /* solution will be stored here */
166c4762a1bSJed Brown   Mat                A; /* Jacobian matrix */
167c4762a1bSJed Brown   PetscInt           n = 2;
168c4762a1bSJed Brown   PetscReal          ftime;
169c4762a1bSJed Brown   PetscInt           steps;
170c4762a1bSJed Brown   PetscScalar       *u;
171c4762a1bSJed Brown   const PetscScalar *x_ptr, *qx_ptr;
172c4762a1bSJed Brown   Vec                q;
173c4762a1bSJed Brown   PetscInt           direction[2];
174c4762a1bSJed Brown   PetscBool          terminate[2];
175c4762a1bSJed Brown 
176*3ba16761SJacob Faibussowitsch   PetscFunctionBeginUser;
1779566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(P, &x_ptr));
178c4762a1bSJed Brown   ctx->Pm = x_ptr[0];
1799566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(P, &x_ptr));
180c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
181c4762a1bSJed Brown     Create necessary matrix and vectors
182c4762a1bSJed Brown     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1839566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
1849566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
1859566063dSJacob Faibussowitsch   PetscCall(MatSetType(A, MATDENSE));
1869566063dSJacob Faibussowitsch   PetscCall(MatSetFromOptions(A));
1879566063dSJacob Faibussowitsch   PetscCall(MatSetUp(A));
188c4762a1bSJed Brown 
1899566063dSJacob Faibussowitsch   PetscCall(MatCreateVecs(A, &U, NULL));
190c4762a1bSJed Brown 
191c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
192c4762a1bSJed Brown      Create timestepping solver context
193c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
1949566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1959566063dSJacob Faibussowitsch   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
1969566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSCN));
1979566063dSJacob Faibussowitsch   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, ctx));
1989566063dSJacob Faibussowitsch   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, ctx));
199c4762a1bSJed Brown 
200c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
201c4762a1bSJed Brown      Set initial conditions
202c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2039566063dSJacob Faibussowitsch   PetscCall(VecGetArray(U, &u));
204c4762a1bSJed Brown   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
205c4762a1bSJed Brown   u[1] = 1.0;
2069566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(U, &u));
2079566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, U));
208c4762a1bSJed Brown 
209c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
210c4762a1bSJed Brown      Set solver options
211c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2129566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, 1.0));
2139566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
2149566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.03125));
2159566063dSJacob Faibussowitsch   PetscCall(TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts));
2169566063dSJacob Faibussowitsch   PetscCall(TSGetSolution(quadts, &q));
2179566063dSJacob Faibussowitsch   PetscCall(VecSet(q, 0.0));
2189566063dSJacob Faibussowitsch   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, ctx));
2199566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
220c4762a1bSJed Brown 
221c4762a1bSJed Brown   direction[0] = direction[1] = 1;
222c4762a1bSJed Brown   terminate[0] = terminate[1] = PETSC_FALSE;
223c4762a1bSJed Brown 
2249566063dSJacob Faibussowitsch   PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)ctx));
225c4762a1bSJed Brown 
226c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
227c4762a1bSJed Brown      Solve nonlinear system
228c4762a1bSJed Brown      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2299566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, U));
230c4762a1bSJed Brown 
2319566063dSJacob Faibussowitsch   PetscCall(TSGetSolveTime(ts, &ftime));
2329566063dSJacob Faibussowitsch   PetscCall(TSGetStepNumber(ts, &steps));
2339566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(q, &qx_ptr));
234c4762a1bSJed Brown   *f = -ctx->Pm + qx_ptr[0];
2359566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(q, &qx_ptr));
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
238c4762a1bSJed Brown      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
239c4762a1bSJed Brown    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
2409566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&A));
2419566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&U));
2429566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
243*3ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
244c4762a1bSJed Brown }
245c4762a1bSJed Brown 
246c4762a1bSJed Brown /*TEST
247c4762a1bSJed Brown 
248c4762a1bSJed Brown    build:
249c4762a1bSJed Brown       requires: !complex !single
250c4762a1bSJed Brown 
251c4762a1bSJed Brown    test:
252c4762a1bSJed Brown       args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3
253c4762a1bSJed Brown 
254c4762a1bSJed Brown TEST*/
255