Actual source code: ex3opt_fd.c

  1: static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\n";

  3: /*F

  5: \begin{eqnarray}
  6:                  \frac{d \theta}{dt} = \omega_b (\omega - \omega_s)
  7:                  \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\
  8: \end{eqnarray}

 10: F*/

 12: /*
 13:   Solve the same optimization problem as in ex3opt.c.
 14:   Use finite difference to approximate the gradients.
 15: */
 16: #include <petsctao.h>
 17: #include <petscts.h>
 18: #include "ex3.h"

 20: PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *);

 22: PetscErrorCode monitor(Tao tao, AppCtx *ctx)
 23: {
 24:   FILE              *fp;
 25:   PetscInt           iterate;
 26:   PetscReal          f, gnorm, cnorm, xdiff;
 27:   Vec                X, G;
 28:   const PetscScalar *x, *g;
 29:   TaoConvergedReason reason;

 31:   PetscFunctionBeginUser;
 32:   PetscCall(TaoGetSolutionStatus(tao, &iterate, &f, &gnorm, &cnorm, &xdiff, &reason));
 33:   PetscCall(TaoGetSolution(tao, &X));
 34:   PetscCall(TaoGetGradient(tao, &G, NULL, NULL));
 35:   PetscCall(VecGetArrayRead(X, &x));
 36:   PetscCall(VecGetArrayRead(G, &g));
 37:   fp = fopen("ex3opt_fd_conv.out", "a");
 38:   PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fp, "%" PetscInt_FMT " %g %.12lf %.12lf\n", iterate, (double)gnorm, (double)PetscRealPart(x[0]), (double)PetscRealPart(g[0])));
 39:   PetscCall(VecRestoreArrayRead(X, &x));
 40:   PetscCall(VecRestoreArrayRead(G, &g));
 41:   fclose(fp);
 42:   PetscFunctionReturn(PETSC_SUCCESS);
 43: }

 45: int main(int argc, char **argv)
 46: {
 47:   Vec          p;
 48:   PetscScalar *x_ptr;
 49:   PetscMPIInt  size;
 50:   AppCtx       ctx;
 51:   Vec          lowerb, upperb;
 52:   Tao          tao;
 53:   KSP          ksp;
 54:   PC           pc;
 55:   PetscBool    printtofile;
 56:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 57:      Initialize program
 58:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 59:   PetscFunctionBeginUser;
 60:   PetscCall(PetscInitialize(&argc, &argv, NULL, help));
 61:   PetscFunctionBeginUser;
 62:   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
 63:   PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!");

 65:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
 66:     Set runtime options
 67:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
 68:   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", "");
 69:   {
 70:     ctx.beta    = 2;
 71:     ctx.c       = 10000.0;
 72:     ctx.u_s     = 1.0;
 73:     ctx.omega_s = 1.0;
 74:     ctx.omega_b = 120.0 * PETSC_PI;
 75:     ctx.H       = 5.0;
 76:     PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL));
 77:     ctx.D = 5.0;
 78:     PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL));
 79:     ctx.E        = 1.1378;
 80:     ctx.V        = 1.0;
 81:     ctx.X        = 0.545;
 82:     ctx.Pmax     = ctx.E * ctx.V / ctx.X;
 83:     ctx.Pmax_ini = ctx.Pmax;
 84:     PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL));
 85:     ctx.Pm = 1.06;
 86:     PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL));
 87:     ctx.tf  = 0.1;
 88:     ctx.tcl = 0.2;
 89:     PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL));
 90:     PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL));
 91:     printtofile = PETSC_FALSE;
 92:     PetscCall(PetscOptionsBool("-printtofile", "Print convergence results to file", "", printtofile, &printtofile, NULL));
 93:   }
 94:   PetscOptionsEnd();

 96:   /* Create TAO solver and set desired solution method */
 97:   PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao));
 98:   PetscCall(TaoSetType(tao, TAOBLMVM));
 99:   if (printtofile) PetscCall(TaoSetMonitor(tao, (PetscErrorCode(*)(Tao, void *))monitor, (void *)&ctx, NULL));
100:   PetscCall(TaoSetMaximumIterations(tao, 30));
101:   /*
102:      Optimization starts
103:   */
104:   /* Set initial solution guess */
105:   PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p));
106:   PetscCall(VecGetArray(p, &x_ptr));
107:   x_ptr[0] = ctx.Pm;
108:   PetscCall(VecRestoreArray(p, &x_ptr));

110:   PetscCall(TaoSetSolution(tao, p));
111:   /* Set routine for function and gradient evaluation */
112:   PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx));
113:   PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, (void *)&ctx));

115:   /* Set bounds for the optimization */
116:   PetscCall(VecDuplicate(p, &lowerb));
117:   PetscCall(VecDuplicate(p, &upperb));
118:   PetscCall(VecGetArray(lowerb, &x_ptr));
119:   x_ptr[0] = 0.;
120:   PetscCall(VecRestoreArray(lowerb, &x_ptr));
121:   PetscCall(VecGetArray(upperb, &x_ptr));
122:   x_ptr[0] = 1.1;
123:   PetscCall(VecRestoreArray(upperb, &x_ptr));
124:   PetscCall(TaoSetVariableBounds(tao, lowerb, upperb));

126:   /* Check for any TAO command line options */
127:   PetscCall(TaoSetFromOptions(tao));
128:   PetscCall(TaoGetKSP(tao, &ksp));
129:   if (ksp) {
130:     PetscCall(KSPGetPC(ksp, &pc));
131:     PetscCall(PCSetType(pc, PCNONE));
132:   }

134:   /* SOLVE THE APPLICATION */
135:   PetscCall(TaoSolve(tao));

137:   PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD));

139:   /* Free TAO data structures */
140:   PetscCall(TaoDestroy(&tao));
141:   PetscCall(VecDestroy(&p));
142:   PetscCall(VecDestroy(&lowerb));
143:   PetscCall(VecDestroy(&upperb));
144:   PetscCall(PetscFinalize());
145:   return 0;
146: }

148: /* ------------------------------------------------------------------ */
149: /*
150:    FormFunction - Evaluates the function and corresponding gradient.

152:    Input Parameters:
153:    tao - the Tao context
154:    X   - the input vector
155:    ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient()

157:    Output Parameters:
158:    f   - the newly evaluated function
159: */
160: PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0)
161: {
162:   AppCtx            *ctx = (AppCtx *)ctx0;
163:   TS                 ts, quadts;
164:   Vec                U; /* solution will be stored here */
165:   Mat                A; /* Jacobian matrix */
166:   PetscInt           n = 2;
167:   PetscReal          ftime;
168:   PetscInt           steps;
169:   PetscScalar       *u;
170:   const PetscScalar *x_ptr, *qx_ptr;
171:   Vec                q;
172:   PetscInt           direction[2];
173:   PetscBool          terminate[2];

175:   PetscFunctionBeginUser;
176:   PetscCall(VecGetArrayRead(P, &x_ptr));
177:   ctx->Pm = x_ptr[0];
178:   PetscCall(VecRestoreArrayRead(P, &x_ptr));
179:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
180:     Create necessary matrix and vectors
181:     - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
182:   PetscCall(MatCreate(PETSC_COMM_WORLD, &A));
183:   PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE));
184:   PetscCall(MatSetType(A, MATDENSE));
185:   PetscCall(MatSetFromOptions(A));
186:   PetscCall(MatSetUp(A));

188:   PetscCall(MatCreateVecs(A, &U, NULL));

190:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
191:      Create timestepping solver context
192:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
193:   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
194:   PetscCall(TSSetProblemType(ts, TS_NONLINEAR));
195:   PetscCall(TSSetType(ts, TSCN));
196:   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, ctx));
197:   PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, ctx));

199:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
200:      Set initial conditions
201:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
202:   PetscCall(VecGetArray(U, &u));
203:   u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax);
204:   u[1] = 1.0;
205:   PetscCall(VecRestoreArray(U, &u));
206:   PetscCall(TSSetSolution(ts, U));

208:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
209:      Set solver options
210:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
211:   PetscCall(TSSetMaxTime(ts, 1.0));
212:   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP));
213:   PetscCall(TSSetTimeStep(ts, 0.03125));
214:   PetscCall(TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts));
215:   PetscCall(TSGetSolution(quadts, &q));
216:   PetscCall(VecSet(q, 0.0));
217:   PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, ctx));
218:   PetscCall(TSSetFromOptions(ts));

220:   direction[0] = direction[1] = 1;
221:   terminate[0] = terminate[1] = PETSC_FALSE;

223:   PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)ctx));

225:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
226:      Solve nonlinear system
227:      - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
228:   PetscCall(TSSolve(ts, U));

230:   PetscCall(TSGetSolveTime(ts, &ftime));
231:   PetscCall(TSGetStepNumber(ts, &steps));
232:   PetscCall(VecGetArrayRead(q, &qx_ptr));
233:   *f = -ctx->Pm + qx_ptr[0];
234:   PetscCall(VecRestoreArrayRead(q, &qx_ptr));

236:   /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - -
237:      Free work space.  All PETSc objects should be destroyed when they are no longer needed.
238:    - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */
239:   PetscCall(MatDestroy(&A));
240:   PetscCall(VecDestroy(&U));
241:   PetscCall(TSDestroy(&ts));
242:   PetscFunctionReturn(PETSC_SUCCESS);
243: }

245: /*TEST

247:    build:
248:       requires: !complex !single

250:    test:
251:       args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3

253: TEST*/