1c4762a1bSJed Brown static char help[] = "Finds optimal parameter P_m for the generator system while maintaining generator stability.\n"; 2c4762a1bSJed Brown 3c4762a1bSJed Brown /*F 4c4762a1bSJed Brown 5c4762a1bSJed Brown \begin{eqnarray} 6c4762a1bSJed Brown \frac{d \theta}{dt} = \omega_b (\omega - \omega_s) 7c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\omega - \omega_s)\\ 8c4762a1bSJed Brown \end{eqnarray} 9c4762a1bSJed Brown 10c4762a1bSJed Brown F*/ 11c4762a1bSJed Brown 12c4762a1bSJed Brown /* 13c4762a1bSJed Brown Solve the same optimization problem as in ex3opt.c. 14c4762a1bSJed Brown Use finite difference to approximate the gradients. 15c4762a1bSJed Brown */ 16c4762a1bSJed Brown #include <petsctao.h> 17c4762a1bSJed Brown #include <petscts.h> 18c4762a1bSJed Brown #include "ex3.h" 19c4762a1bSJed Brown 20c4762a1bSJed Brown PetscErrorCode FormFunction(Tao, Vec, PetscReal *, void *); 21c4762a1bSJed Brown 22d71ae5a4SJacob Faibussowitsch PetscErrorCode monitor(Tao tao, AppCtx *ctx) 23d71ae5a4SJacob Faibussowitsch { 24c4762a1bSJed Brown FILE *fp; 25c4762a1bSJed Brown PetscInt iterate; 26c4762a1bSJed Brown PetscReal f, gnorm, cnorm, xdiff; 27c4762a1bSJed Brown Vec X, G; 28c4762a1bSJed Brown const PetscScalar *x, *g; 29c4762a1bSJed Brown TaoConvergedReason reason; 30c4762a1bSJed Brown 31c4762a1bSJed Brown PetscFunctionBeginUser; 329566063dSJacob Faibussowitsch PetscCall(TaoGetSolutionStatus(tao, &iterate, &f, &gnorm, &cnorm, &xdiff, &reason)); 339566063dSJacob Faibussowitsch PetscCall(TaoGetSolution(tao, &X)); 349566063dSJacob Faibussowitsch PetscCall(TaoGetGradient(tao, &G, NULL, NULL)); 359566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(X, &x)); 369566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(G, &g)); 37c4762a1bSJed Brown fp = fopen("ex3opt_fd_conv.out", "a"); 3863a3b9bcSJacob Faibussowitsch PetscCall(PetscFPrintf(PETSC_COMM_WORLD, fp, "%" PetscInt_FMT " %g %.12lf %.12lf\n", iterate, (double)gnorm, (double)PetscRealPart(x[0]), (double)PetscRealPart(g[0]))); 399566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, &x)); 409566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(G, &g)); 41c4762a1bSJed Brown fclose(fp); 423ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 43c4762a1bSJed Brown } 44c4762a1bSJed Brown 45d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 46d71ae5a4SJacob Faibussowitsch { 47c4762a1bSJed Brown Vec p; 48c4762a1bSJed Brown PetscScalar *x_ptr; 49c4762a1bSJed Brown PetscMPIInt size; 50c4762a1bSJed Brown AppCtx ctx; 51c4762a1bSJed Brown Vec lowerb, upperb; 52c4762a1bSJed Brown Tao tao; 53c4762a1bSJed Brown KSP ksp; 54c4762a1bSJed Brown PC pc; 55c4762a1bSJed Brown PetscBool printtofile; 56c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 57c4762a1bSJed Brown Initialize program 58c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 59327415f7SBarry Smith PetscFunctionBeginUser; 609566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, NULL, help)); 619566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 623c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "This is a uniprocessor example only!"); 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 65c4762a1bSJed Brown Set runtime options 66c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 67d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 68c4762a1bSJed Brown { 69c4762a1bSJed Brown ctx.beta = 2; 70c4762a1bSJed Brown ctx.c = 10000.0; 71c4762a1bSJed Brown ctx.u_s = 1.0; 72c4762a1bSJed Brown ctx.omega_s = 1.0; 73c4762a1bSJed Brown ctx.omega_b = 120.0 * PETSC_PI; 74c4762a1bSJed Brown ctx.H = 5.0; 759566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 76c4762a1bSJed Brown ctx.D = 5.0; 779566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 78c4762a1bSJed Brown ctx.E = 1.1378; 79c4762a1bSJed Brown ctx.V = 1.0; 80c4762a1bSJed Brown ctx.X = 0.545; 81c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 82c4762a1bSJed Brown ctx.Pmax_ini = ctx.Pmax; 839566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 84c4762a1bSJed Brown ctx.Pm = 1.06; 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 86c4762a1bSJed Brown ctx.tf = 0.1; 87c4762a1bSJed Brown ctx.tcl = 0.2; 889566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 899566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 90c4762a1bSJed Brown printtofile = PETSC_FALSE; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-printtofile", "Print convergence results to file", "", printtofile, &printtofile, NULL)); 92c4762a1bSJed Brown } 93d0609cedSBarry Smith PetscOptionsEnd(); 94c4762a1bSJed Brown 95c4762a1bSJed Brown /* Create TAO solver and set desired solution method */ 969566063dSJacob Faibussowitsch PetscCall(TaoCreate(PETSC_COMM_WORLD, &tao)); 979566063dSJacob Faibussowitsch PetscCall(TaoSetType(tao, TAOBLMVM)); 9810978b7dSBarry Smith if (printtofile) PetscCall(TaoMonitorSet(tao, (PetscErrorCode (*)(Tao, void *))monitor, (void *)&ctx, NULL)); 999566063dSJacob Faibussowitsch PetscCall(TaoSetMaximumIterations(tao, 30)); 100c4762a1bSJed Brown /* 101c4762a1bSJed Brown Optimization starts 102c4762a1bSJed Brown */ 103c4762a1bSJed Brown /* Set initial solution guess */ 1049566063dSJacob Faibussowitsch PetscCall(VecCreateSeq(PETSC_COMM_WORLD, 1, &p)); 1059566063dSJacob Faibussowitsch PetscCall(VecGetArray(p, &x_ptr)); 106c4762a1bSJed Brown x_ptr[0] = ctx.Pm; 1079566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(p, &x_ptr)); 108c4762a1bSJed Brown 1099566063dSJacob Faibussowitsch PetscCall(TaoSetSolution(tao, p)); 110c4762a1bSJed Brown /* Set routine for function and gradient evaluation */ 1119566063dSJacob Faibussowitsch PetscCall(TaoSetObjective(tao, FormFunction, (void *)&ctx)); 1129566063dSJacob Faibussowitsch PetscCall(TaoSetGradient(tao, NULL, TaoDefaultComputeGradient, (void *)&ctx)); 113c4762a1bSJed Brown 114c4762a1bSJed Brown /* Set bounds for the optimization */ 1159566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &lowerb)); 1169566063dSJacob Faibussowitsch PetscCall(VecDuplicate(p, &upperb)); 1179566063dSJacob Faibussowitsch PetscCall(VecGetArray(lowerb, &x_ptr)); 118c4762a1bSJed Brown x_ptr[0] = 0.; 1199566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(lowerb, &x_ptr)); 1209566063dSJacob Faibussowitsch PetscCall(VecGetArray(upperb, &x_ptr)); 121c4762a1bSJed Brown x_ptr[0] = 1.1; 1229566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(upperb, &x_ptr)); 1239566063dSJacob Faibussowitsch PetscCall(TaoSetVariableBounds(tao, lowerb, upperb)); 124c4762a1bSJed Brown 125c4762a1bSJed Brown /* Check for any TAO command line options */ 1269566063dSJacob Faibussowitsch PetscCall(TaoSetFromOptions(tao)); 1279566063dSJacob Faibussowitsch PetscCall(TaoGetKSP(tao, &ksp)); 128c4762a1bSJed Brown if (ksp) { 1299566063dSJacob Faibussowitsch PetscCall(KSPGetPC(ksp, &pc)); 1309566063dSJacob Faibussowitsch PetscCall(PCSetType(pc, PCNONE)); 131c4762a1bSJed Brown } 132c4762a1bSJed Brown 133c4762a1bSJed Brown /* SOLVE THE APPLICATION */ 1349566063dSJacob Faibussowitsch PetscCall(TaoSolve(tao)); 135c4762a1bSJed Brown 1369566063dSJacob Faibussowitsch PetscCall(VecView(p, PETSC_VIEWER_STDOUT_WORLD)); 137c4762a1bSJed Brown 138c4762a1bSJed Brown /* Free TAO data structures */ 1399566063dSJacob Faibussowitsch PetscCall(TaoDestroy(&tao)); 1409566063dSJacob Faibussowitsch PetscCall(VecDestroy(&p)); 1419566063dSJacob Faibussowitsch PetscCall(VecDestroy(&lowerb)); 1429566063dSJacob Faibussowitsch PetscCall(VecDestroy(&upperb)); 1439566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 144b122ec5aSJacob Faibussowitsch return 0; 145c4762a1bSJed Brown } 146c4762a1bSJed Brown 147c4762a1bSJed Brown /* ------------------------------------------------------------------ */ 148c4762a1bSJed Brown /* 149c4762a1bSJed Brown FormFunction - Evaluates the function and corresponding gradient. 150c4762a1bSJed Brown 151c4762a1bSJed Brown Input Parameters: 152c4762a1bSJed Brown tao - the Tao context 153c4762a1bSJed Brown X - the input vector 154a82e8c82SStefano Zampini ptr - optional user-defined context, as set by TaoSetObjectiveAndGradient() 155c4762a1bSJed Brown 156c4762a1bSJed Brown Output Parameters: 157c4762a1bSJed Brown f - the newly evaluated function 158c4762a1bSJed Brown */ 159*2a8381b2SBarry Smith PetscErrorCode FormFunction(Tao tao, Vec P, PetscReal *f, PetscCtx ctx0) 160d71ae5a4SJacob Faibussowitsch { 161c4762a1bSJed Brown AppCtx *ctx = (AppCtx *)ctx0; 162c4762a1bSJed Brown TS ts, quadts; 163c4762a1bSJed Brown Vec U; /* solution will be stored here */ 164c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 165c4762a1bSJed Brown PetscInt n = 2; 166c4762a1bSJed Brown PetscReal ftime; 167c4762a1bSJed Brown PetscInt steps; 168c4762a1bSJed Brown PetscScalar *u; 169c4762a1bSJed Brown const PetscScalar *x_ptr, *qx_ptr; 170c4762a1bSJed Brown Vec q; 171c4762a1bSJed Brown PetscInt direction[2]; 172c4762a1bSJed Brown PetscBool terminate[2]; 173c4762a1bSJed Brown 1743ba16761SJacob Faibussowitsch PetscFunctionBeginUser; 1759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(P, &x_ptr)); 176c4762a1bSJed Brown ctx->Pm = x_ptr[0]; 1779566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(P, &x_ptr)); 178c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 179c4762a1bSJed Brown Create necessary matrix and vectors 180c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1819566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1829566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1839566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 1849566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1859566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 186c4762a1bSJed Brown 1879566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 190c4762a1bSJed Brown Create timestepping solver context 191c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1929566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1939566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1949566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSCN)); 1958434afd1SBarry Smith PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, ctx)); 1968434afd1SBarry Smith PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, ctx)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 199c4762a1bSJed Brown Set initial conditions 200c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2019566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 202c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx->Pm / ctx->Pmax); 203c4762a1bSJed Brown u[1] = 1.0; 2049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 2059566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 208c4762a1bSJed Brown Set solver options 209c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2109566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 1.0)); 2119566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 2129566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, 0.03125)); 2139566063dSJacob Faibussowitsch PetscCall(TSCreateQuadratureTS(ts, PETSC_TRUE, &quadts)); 2149566063dSJacob Faibussowitsch PetscCall(TSGetSolution(quadts, &q)); 2159566063dSJacob Faibussowitsch PetscCall(VecSet(q, 0.0)); 2168434afd1SBarry Smith PetscCall(TSSetRHSFunction(quadts, NULL, (TSRHSFunctionFn *)CostIntegrand, ctx)); 2179566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 218c4762a1bSJed Brown 219c4762a1bSJed Brown direction[0] = direction[1] = 1; 220c4762a1bSJed Brown terminate[0] = terminate[1] = PETSC_FALSE; 221c4762a1bSJed Brown 2229566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)ctx)); 223c4762a1bSJed Brown 224c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225c4762a1bSJed Brown Solve nonlinear system 226c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2279566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 228c4762a1bSJed Brown 2299566063dSJacob Faibussowitsch PetscCall(TSGetSolveTime(ts, &ftime)); 2309566063dSJacob Faibussowitsch PetscCall(TSGetStepNumber(ts, &steps)); 2319566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(q, &qx_ptr)); 232c4762a1bSJed Brown *f = -ctx->Pm + qx_ptr[0]; 2339566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(q, &qx_ptr)); 234c4762a1bSJed Brown 235c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 236c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 237c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2389566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2399566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 2409566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2413ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 242c4762a1bSJed Brown } 243c4762a1bSJed Brown 244c4762a1bSJed Brown /*TEST 245c4762a1bSJed Brown 246c4762a1bSJed Brown build: 247c4762a1bSJed Brown requires: !complex !single 248c4762a1bSJed Brown 249c4762a1bSJed Brown test: 250c4762a1bSJed Brown args: -ts_type cn -pc_type lu -tao_monitor -tao_gatol 1e-3 251c4762a1bSJed Brown 252c4762a1bSJed Brown TEST*/ 253