1c4762a1bSJed Brown 2c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\n"; 3c4762a1bSJed Brown 4c4762a1bSJed Brown /*F 5c4762a1bSJed Brown 6c4762a1bSJed Brown \begin{eqnarray} 7c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - \frac{EV}{X} \sin(\theta) -D(\omega - \omega_s)\\ 8c4762a1bSJed Brown \frac{d \theta}{dt} = \omega - \omega_s 9c4762a1bSJed Brown \end{eqnarray} 10c4762a1bSJed Brown 11c4762a1bSJed Brown Ensemble of initial conditions 12c4762a1bSJed Brown ./ex2 -ensemble -ts_monitor_draw_solution_phase -1,-3,3,3 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 13c4762a1bSJed Brown 14c4762a1bSJed Brown Fault at .1 seconds 15c4762a1bSJed Brown ./ex2 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 16c4762a1bSJed Brown 17c4762a1bSJed Brown Initial conditions same as when fault is ended 18c4762a1bSJed Brown ./ex2 -u 0.496792,1.00932 -ts_monitor_draw_solution_phase .42,.95,.6,1.05 -ts_adapt_dt_max .01 -ts_monitor -ts_type rosw -pc_type lu -ksp_type preonly 19c4762a1bSJed Brown 20c4762a1bSJed Brown F*/ 21c4762a1bSJed Brown 22c4762a1bSJed Brown /* 23c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 24c4762a1bSJed Brown file automatically includes: 25c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 26c4762a1bSJed Brown petscmat.h - matrices 27c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 28c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 29c4762a1bSJed Brown petscksp.h - linear solvers 30c4762a1bSJed Brown */ 31c4762a1bSJed Brown 32c4762a1bSJed Brown #include <petscts.h> 33c4762a1bSJed Brown 34c4762a1bSJed Brown typedef struct { 35c4762a1bSJed Brown PetscScalar H, D, omega_s, Pmax, Pm, E, V, X; 36c4762a1bSJed Brown PetscReal tf, tcl; 37c4762a1bSJed Brown } AppCtx; 38c4762a1bSJed Brown 39c4762a1bSJed Brown /* 40c4762a1bSJed Brown Defines the ODE passed to the ODE solver 41c4762a1bSJed Brown */ 42*9371c9d4SSatish Balay static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) { 43c4762a1bSJed Brown PetscScalar *f, Pmax; 44c4762a1bSJed Brown const PetscScalar *u, *udot; 45c4762a1bSJed Brown 46c4762a1bSJed Brown PetscFunctionBegin; 47c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 489566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 509566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 51c4762a1bSJed Brown if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */ 52c4762a1bSJed Brown else if (t >= ctx->tcl) Pmax = ctx->E / 0.745; 53c4762a1bSJed Brown else Pmax = ctx->Pmax; 54c4762a1bSJed Brown f[0] = udot[0] - ctx->omega_s * (u[1] - 1.0); 55c4762a1bSJed Brown f[1] = 2.0 * ctx->H * udot[1] + Pmax * PetscSinScalar(u[0]) + ctx->D * (u[1] - 1.0) - ctx->Pm; 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 599566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 60c4762a1bSJed Brown PetscFunctionReturn(0); 61c4762a1bSJed Brown } 62c4762a1bSJed Brown 63c4762a1bSJed Brown /* 64c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 65c4762a1bSJed Brown */ 66*9371c9d4SSatish Balay static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx) { 67c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 68c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 69c4762a1bSJed Brown const PetscScalar *u, *udot; 70c4762a1bSJed Brown 71c4762a1bSJed Brown PetscFunctionBegin; 729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 739566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 74c4762a1bSJed Brown if ((t > ctx->tf) && (t < ctx->tcl)) Pmax = 0.0; /* A short-circuit on the generator terminal that drives the electrical power output (Pmax*sin(delta)) to 0 */ 75c4762a1bSJed Brown else if (t >= ctx->tcl) Pmax = ctx->E / 0.745; 76c4762a1bSJed Brown else Pmax = ctx->Pmax; 77c4762a1bSJed Brown 78*9371c9d4SSatish Balay J[0][0] = a; 79*9371c9d4SSatish Balay J[0][1] = -ctx->omega_s; 80*9371c9d4SSatish Balay J[1][1] = 2.0 * ctx->H * a + ctx->D; 81*9371c9d4SSatish Balay J[1][0] = Pmax * PetscCosScalar(u[0]); 82c4762a1bSJed Brown 839566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 849566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 859566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 86c4762a1bSJed Brown 879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 89c4762a1bSJed Brown if (A != B) { 909566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 919566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 92c4762a1bSJed Brown } 93c4762a1bSJed Brown PetscFunctionReturn(0); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 96*9371c9d4SSatish Balay PetscErrorCode PostStep(TS ts) { 97c4762a1bSJed Brown Vec X; 98c4762a1bSJed Brown PetscReal t; 99c4762a1bSJed Brown 100c4762a1bSJed Brown PetscFunctionBegin; 1019566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 102c4762a1bSJed Brown if (t >= .2) { 1039566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 1049566063dSJacob Faibussowitsch PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 105c4762a1bSJed Brown exit(0); 106c4762a1bSJed Brown /* results in initial conditions after fault of -u 0.496792,1.00932 */ 107c4762a1bSJed Brown } 108c4762a1bSJed Brown PetscFunctionReturn(0); 109c4762a1bSJed Brown } 110c4762a1bSJed Brown 111*9371c9d4SSatish Balay int main(int argc, char **argv) { 112c4762a1bSJed Brown TS ts; /* ODE integrator */ 113c4762a1bSJed Brown Vec U; /* solution will be stored here */ 114c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 115c4762a1bSJed Brown PetscMPIInt size; 116c4762a1bSJed Brown PetscInt n = 2; 117c4762a1bSJed Brown AppCtx ctx; 118c4762a1bSJed Brown PetscScalar *u; 119c4762a1bSJed Brown PetscReal du[2] = {0.0, 0.0}; 120c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE, flg1, flg2; 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123c4762a1bSJed Brown Initialize program 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125327415f7SBarry Smith PetscFunctionBeginUser; 1269566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1279566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1283c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 129c4762a1bSJed Brown 130c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 131c4762a1bSJed Brown Create necessary matrix and vectors 132c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1339566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1349566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1359566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 1369566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1379566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142c4762a1bSJed Brown Set runtime options 143c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 145c4762a1bSJed Brown { 146c4762a1bSJed Brown ctx.omega_s = 2.0 * PETSC_PI * 60.0; 147c4762a1bSJed Brown ctx.H = 5.0; 1489566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 149c4762a1bSJed Brown ctx.D = 5.0; 1509566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 151c4762a1bSJed Brown ctx.E = 1.1378; 152c4762a1bSJed Brown ctx.V = 1.0; 153c4762a1bSJed Brown ctx.X = 0.545; 154c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 1559566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 156c4762a1bSJed Brown ctx.Pm = 0.9; 1579566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 158c4762a1bSJed Brown ctx.tf = 1.0; 159c4762a1bSJed Brown ctx.tcl = 1.05; 1609566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 1619566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 1629566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL)); 163c4762a1bSJed Brown if (ensemble) { 164c4762a1bSJed Brown ctx.tf = -1; 165c4762a1bSJed Brown ctx.tcl = -1; 166c4762a1bSJed Brown } 167c4762a1bSJed Brown 1689566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 169c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 170c4762a1bSJed Brown u[1] = 1.0; 1719566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1)); 172c4762a1bSJed Brown n = 2; 1739566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2)); 174c4762a1bSJed Brown u[0] += du[0]; 175c4762a1bSJed Brown u[1] += du[1]; 1769566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 177c4762a1bSJed Brown if (flg1 || flg2) { 178c4762a1bSJed Brown ctx.tf = -1; 179c4762a1bSJed Brown ctx.tcl = -1; 180c4762a1bSJed Brown } 181c4762a1bSJed Brown } 182d0609cedSBarry Smith PetscOptionsEnd(); 183c4762a1bSJed Brown 184c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 185c4762a1bSJed Brown Create timestepping solver context 186c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1879566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1889566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1899566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSROSW)); 1909566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &ctx)); 1919566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx)); 192c4762a1bSJed Brown 193c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 194c4762a1bSJed Brown Set initial conditions 195c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1969566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 199c4762a1bSJed Brown Set solver options 200c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2019566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 35.0)); 2029566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 2039566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 2049566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 2059566063dSJacob Faibussowitsch /* PetscCall(TSSetPostStep(ts,PostStep)); */ 206c4762a1bSJed Brown 207c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 208c4762a1bSJed Brown Solve nonlinear system 209c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 210c4762a1bSJed Brown if (ensemble) { 211c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 2129566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 213c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 214c4762a1bSJed Brown u[1] = ctx.omega_s; 215c4762a1bSJed Brown u[0] += du[0]; 216c4762a1bSJed Brown u[1] += du[1]; 2179566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 2189566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 2199566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 220c4762a1bSJed Brown } 221c4762a1bSJed Brown } else { 2229566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 223c4762a1bSJed Brown } 224c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 225c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 226c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2279566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2289566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 2299566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2309566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 231b122ec5aSJacob Faibussowitsch return 0; 232c4762a1bSJed Brown } 233c4762a1bSJed Brown 234c4762a1bSJed Brown /*TEST 235c4762a1bSJed Brown 236c4762a1bSJed Brown build: 237c4762a1bSJed Brown requires: !complex 238c4762a1bSJed Brown 239c4762a1bSJed Brown test: 240c4762a1bSJed Brown args: -nox -ts_dt 10 241c4762a1bSJed Brown 242c4762a1bSJed Brown TEST*/ 243