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 FILE *fp; 25 PetscInt iterate; 26 PetscReal f, gnorm, cnorm, xdiff; 27 Vec X, G; 28 const PetscScalar *x, *g; 29 TaoConvergedReason reason; 30 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(0); 43 } 44 45 int main(int argc, char **argv) { 46 Vec p; 47 PetscScalar *x_ptr; 48 PetscMPIInt size; 49 AppCtx ctx; 50 Vec lowerb, upperb; 51 Tao tao; 52 KSP ksp; 53 PC pc; 54 PetscBool printtofile; 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 } 93 PetscOptionsEnd(); 94 95 /* Create TAO solver and set desired solution method */ 96 PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 97 PetscCall(TaoSetType(tao, TAOBLMVM)); 98 if (printtofile) { PetscCall(TaoSetMonitor(tao, (PetscErrorCode(*)(Tao, void *))monitor, (void *)&ctx, PETSC_NULL)); } 99 PetscCall(TaoSetMaximumIterations(tao, 30)); 100 /* 101 Optimization starts 102 */ 103 /* Set initial solution guess */ 104 PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p)); 105 PetscCall(VecGetArray(p, &x_ptr)); 106 x_ptr[0] = ctx.Pm; 107 PetscCall(VecRestoreArray(p, &x_ptr)); 108 109 PetscCall(TaoSetSolution(tao, p)); 110 /* Set routine for function and gradient evaluation */ 111 PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx)); 112 PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, (void *)&ctx)); 113 114 /* Set bounds for the optimization */ 115 PetscCall(VecDuplicate(p, &lowerb)); 116 PetscCall(VecDuplicate(p, &upperb)); 117 PetscCall(VecGetArray(lowerb, &x_ptr)); 118 x_ptr[0] = 0.; 119 PetscCall(VecRestoreArray(lowerb, &x_ptr)); 120 PetscCall(VecGetArray(upperb, &x_ptr)); 121 x_ptr[0] = 1.1; 122 PetscCall(VecRestoreArray(upperb, &x_ptr)); 123 PetscCall(TaoSetVariableBounds(tao, lowerb, upperb)); 124 125 /* Check for any TAO command line options */ 126 PetscCall(TaoSetFromOptions(tao)); 127 PetscCall(TaoGetKSP(tao, &ksp)); 128 if (ksp) { 129 PetscCall(KSPGetPC(ksp, &pc)); 130 PetscCall(PCSetType(pc, PCNONE)); 131 } 132 133 /* SOLVE THE APPLICATION */ 134 PetscCall(TaoSolve(tao)); 135 136 PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD)); 137 138 /* Free TAO data structures */ 139 PetscCall(TaoDestroy(&tao)); 140 PetscCall(VecDestroy(&p)); 141 PetscCall(VecDestroy(&lowerb)); 142 PetscCall(VecDestroy(&upperb)); 143 PetscCall(PetscFinalize()); 144 return 0; 145 } 146 147 /* ------------------------------------------------------------------ */ 148 /* 149 FormFunction - Evaluates the function and corresponding gradient. 150 151 Input Parameters: 152 tao - the Tao context 153 X - the input vector 154 ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 155 156 Output Parameters: 157 f - the newly evaluated function 158 */ 159 PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, void *ctx0) { 160 AppCtx *ctx = (AppCtx *)ctx0; 161 TS ts, quadts; 162 Vec U; /* solution will be stored here */ 163 Mat A; /* Jacobian matrix */ 164 PetscInt n = 2; 165 PetscReal ftime; 166 PetscInt steps; 167 PetscScalar *u; 168 const PetscScalar *x_ptr, *qx_ptr; 169 Vec q; 170 PetscInt direction[2]; 171 PetscBool terminate[2]; 172 173 PetscCall(VecGetArrayRead(P, &x_ptr)); 174 ctx->Pm = x_ptr[0]; 175 PetscCall(VecRestoreArrayRead(P, &x_ptr)); 176 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 177 Create necessary matrix and vectors 178 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 179 PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 180 PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 181 PetscCall(MatSetType(A, MATDENSE)); 182 PetscCall(MatSetFromOptions(A)); 183 PetscCall(MatSetUp(A)); 184 185 PetscCall(MatCreateVecs(A, &U, NULL)); 186 187 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 188 Create timestepping solver context 189 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 190 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 191 PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 192 PetscCall(TSSetType(ts, TSCN)); 193 PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, ctx)); 194 PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, ctx)); 195 196 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 197 Set initial conditions 198 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 199 PetscCall(VecGetArray(U, &u)); 200 u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax); 201 u[1] = 1.0; 202 PetscCall(VecRestoreArray(U, &u)); 203 PetscCall(TSSetSolution(ts, U)); 204 205 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206 Set solver options 207 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 208 PetscCall(TSSetMaxTime(ts, 1.0)); 209 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 210 PetscCall(TSSetTimeStep(ts, 0.03125)); 211 PetscCall(TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts)); 212 PetscCall(TSGetSolution(quadts, &q)); 213 PetscCall(VecSet(q, 0.0)); 214 PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunction)CostIntegrand, ctx)); 215 PetscCall(TSSetFromOptions(ts)); 216 217 direction[0] = direction[1] = 1; 218 terminate[0] = terminate[1] = PETSC_FALSE; 219 220 PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)ctx)); 221 222 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 223 Solve nonlinear system 224 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 225 PetscCall(TSSolve(ts, U)); 226 227 PetscCall(TSGetSolveTime(ts, &ftime)); 228 PetscCall(TSGetStepNumber(ts, &steps)); 229 PetscCall(VecGetArrayRead(q, &qx_ptr)); 230 *f = -ctx->Pm + qx_ptr[0]; 231 PetscCall(VecRestoreArrayRead(q, &qx_ptr)); 232 233 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 234 Free work space. All PETSc objects should be destroyed when they are no longer needed. 235 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 236 PetscCall(MatDestroy(&A)); 237 PetscCall(VecDestroy(&U)); 238 PetscCall(TSDestroy(&ts)); 239 return 0; 240 } 241 242 /*TEST 243 244 build: 245 requires: !complex !single 246 247 test: 248 args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3 249 250 TEST*/ 251