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 */ 42d71ae5a4SJacob Faibussowitsch static PetscErrorCode IFunction(TS ts, PetscReal t, Vec U, Vec Udot, Vec F, AppCtx *ctx) 43d71ae5a4SJacob Faibussowitsch { 44c4762a1bSJed Brown PetscScalar *f, Pmax; 45c4762a1bSJed Brown const PetscScalar *u, *udot; 46c4762a1bSJed Brown 47c4762a1bSJed Brown PetscFunctionBegin; 48c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 499566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 509566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 519566063dSJacob Faibussowitsch PetscCall(VecGetArray(F, &f)); 52c4762a1bSJed 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 */ 53c4762a1bSJed Brown else if (t >= ctx->tcl) Pmax = ctx->E / 0.745; 54c4762a1bSJed Brown else Pmax = ctx->Pmax; 55c4762a1bSJed Brown f[0] = udot[0] - ctx->omega_s * (u[1] - 1.0); 56c4762a1bSJed Brown f[1] = 2.0 * ctx->H * udot[1] + Pmax * PetscSinScalar(u[0]) + ctx->D * (u[1] - 1.0) - ctx->Pm; 57c4762a1bSJed Brown 589566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 599566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 61*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 62c4762a1bSJed Brown } 63c4762a1bSJed Brown 64c4762a1bSJed Brown /* 65c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 66c4762a1bSJed Brown */ 67d71ae5a4SJacob Faibussowitsch static PetscErrorCode IJacobian(TS ts, PetscReal t, Vec U, Vec Udot, PetscReal a, Mat A, Mat B, AppCtx *ctx) 68d71ae5a4SJacob Faibussowitsch { 69c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 70c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 71c4762a1bSJed Brown const PetscScalar *u, *udot; 72c4762a1bSJed Brown 73c4762a1bSJed Brown PetscFunctionBegin; 749566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 759566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(Udot, &udot)); 76c4762a1bSJed 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 */ 77c4762a1bSJed Brown else if (t >= ctx->tcl) Pmax = ctx->E / 0.745; 78c4762a1bSJed Brown else Pmax = ctx->Pmax; 79c4762a1bSJed Brown 809371c9d4SSatish Balay J[0][0] = a; 819371c9d4SSatish Balay J[0][1] = -ctx->omega_s; 829371c9d4SSatish Balay J[1][1] = 2.0 * ctx->H * a + ctx->D; 839371c9d4SSatish Balay J[1][0] = Pmax * PetscCosScalar(u[0]); 84c4762a1bSJed Brown 859566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 869566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 879566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(Udot, &udot)); 88c4762a1bSJed Brown 899566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 909566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 91c4762a1bSJed Brown if (A != B) { 929566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 939566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 94c4762a1bSJed Brown } 95*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98d71ae5a4SJacob Faibussowitsch PetscErrorCode PostStep(TS ts) 99d71ae5a4SJacob Faibussowitsch { 100c4762a1bSJed Brown Vec X; 101c4762a1bSJed Brown PetscReal t; 102c4762a1bSJed Brown 103c4762a1bSJed Brown PetscFunctionBegin; 1049566063dSJacob Faibussowitsch PetscCall(TSGetTime(ts, &t)); 105c4762a1bSJed Brown if (t >= .2) { 1069566063dSJacob Faibussowitsch PetscCall(TSGetSolution(ts, &X)); 1079566063dSJacob Faibussowitsch PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD)); 108c4762a1bSJed Brown exit(0); 109c4762a1bSJed Brown /* results in initial conditions after fault of -u 0.496792,1.00932 */ 110c4762a1bSJed Brown } 111*3ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 112c4762a1bSJed Brown } 113c4762a1bSJed Brown 114d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 115d71ae5a4SJacob Faibussowitsch { 116c4762a1bSJed Brown TS ts; /* ODE integrator */ 117c4762a1bSJed Brown Vec U; /* solution will be stored here */ 118c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 119c4762a1bSJed Brown PetscMPIInt size; 120c4762a1bSJed Brown PetscInt n = 2; 121c4762a1bSJed Brown AppCtx ctx; 122c4762a1bSJed Brown PetscScalar *u; 123c4762a1bSJed Brown PetscReal du[2] = {0.0, 0.0}; 124c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE, flg1, flg2; 125c4762a1bSJed Brown 126c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 127c4762a1bSJed Brown Initialize program 128c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 129327415f7SBarry Smith PetscFunctionBeginUser; 1309566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1319566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1323c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 133c4762a1bSJed Brown 134c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 135c4762a1bSJed Brown Create necessary matrix and vectors 136c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1379566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1389566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1399566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 1409566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1419566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 142c4762a1bSJed Brown 1439566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 144c4762a1bSJed Brown 145c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 146c4762a1bSJed Brown Set runtime options 147c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 148d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 149c4762a1bSJed Brown { 150c4762a1bSJed Brown ctx.omega_s = 2.0 * PETSC_PI * 60.0; 151c4762a1bSJed Brown ctx.H = 5.0; 1529566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 153c4762a1bSJed Brown ctx.D = 5.0; 1549566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 155c4762a1bSJed Brown ctx.E = 1.1378; 156c4762a1bSJed Brown ctx.V = 1.0; 157c4762a1bSJed Brown ctx.X = 0.545; 158c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 1599566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 160c4762a1bSJed Brown ctx.Pm = 0.9; 1619566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 162c4762a1bSJed Brown ctx.tf = 1.0; 163c4762a1bSJed Brown ctx.tcl = 1.05; 1649566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 1659566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 1669566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL)); 167c4762a1bSJed Brown if (ensemble) { 168c4762a1bSJed Brown ctx.tf = -1; 169c4762a1bSJed Brown ctx.tcl = -1; 170c4762a1bSJed Brown } 171c4762a1bSJed Brown 1729566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 173c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 174c4762a1bSJed Brown u[1] = 1.0; 1759566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1)); 176c4762a1bSJed Brown n = 2; 1779566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2)); 178c4762a1bSJed Brown u[0] += du[0]; 179c4762a1bSJed Brown u[1] += du[1]; 1809566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 181c4762a1bSJed Brown if (flg1 || flg2) { 182c4762a1bSJed Brown ctx.tf = -1; 183c4762a1bSJed Brown ctx.tcl = -1; 184c4762a1bSJed Brown } 185c4762a1bSJed Brown } 186d0609cedSBarry Smith PetscOptionsEnd(); 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 189c4762a1bSJed Brown Create timestepping solver context 190c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1919566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1929566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1939566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSROSW)); 1949566063dSJacob Faibussowitsch PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &ctx)); 1959566063dSJacob Faibussowitsch PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobian)IJacobian, &ctx)); 196c4762a1bSJed Brown 197c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 198c4762a1bSJed Brown Set initial conditions 199c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2009566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 201c4762a1bSJed Brown 202c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 203c4762a1bSJed Brown Set solver options 204c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2059566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 35.0)); 2069566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 2079566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 2089566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 2099566063dSJacob Faibussowitsch /* PetscCall(TSSetPostStep(ts,PostStep)); */ 210c4762a1bSJed Brown 211c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 212c4762a1bSJed Brown Solve nonlinear system 213c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 214c4762a1bSJed Brown if (ensemble) { 215c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 2169566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 217c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 218c4762a1bSJed Brown u[1] = ctx.omega_s; 219c4762a1bSJed Brown u[0] += du[0]; 220c4762a1bSJed Brown u[1] += du[1]; 2219566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 2229566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 2239566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 224c4762a1bSJed Brown } 225c4762a1bSJed Brown } else { 2269566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 227c4762a1bSJed Brown } 228c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 229c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 230c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2319566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2329566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 2339566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2349566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 235b122ec5aSJacob Faibussowitsch return 0; 236c4762a1bSJed Brown } 237c4762a1bSJed Brown 238c4762a1bSJed Brown /*TEST 239c4762a1bSJed Brown 240c4762a1bSJed Brown build: 241c4762a1bSJed Brown requires: !complex 242c4762a1bSJed Brown 243c4762a1bSJed Brown test: 244c4762a1bSJed Brown args: -nox -ts_dt 10 245c4762a1bSJed Brown 246c4762a1bSJed Brown TEST*/ 247