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{d \theta}{dt} = \omega_b (\omega - \omega_s) 8c4762a1bSJed Brown \frac{2 H}{\omega_s}\frac{d \omega}{dt} & = & P_m - P_max \sin(\theta) -D(\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_b, 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*d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) 43*d71ae5a4SJacob Faibussowitsch { 44c4762a1bSJed Brown const PetscScalar *u; 45c4762a1bSJed Brown PetscScalar *f, Pmax; 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(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 Pmax = ctx->Pmax; 53c4762a1bSJed Brown 54c4762a1bSJed Brown f[0] = ctx->omega_b * (u[1] - ctx->omega_s); 55c4762a1bSJed Brown f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H); 56c4762a1bSJed Brown 579566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 589566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 59c4762a1bSJed Brown PetscFunctionReturn(0); 60c4762a1bSJed Brown } 61c4762a1bSJed Brown 62c4762a1bSJed Brown /* 63c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 64c4762a1bSJed Brown */ 65*d71ae5a4SJacob Faibussowitsch static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) 66*d71ae5a4SJacob Faibussowitsch { 67c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 68c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 69c4762a1bSJed Brown const PetscScalar *u; 70c4762a1bSJed Brown 71c4762a1bSJed Brown PetscFunctionBegin; 729566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 73c4762a1bSJed 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 */ 74c4762a1bSJed Brown else Pmax = ctx->Pmax; 75c4762a1bSJed Brown 769371c9d4SSatish Balay J[0][0] = 0; 779371c9d4SSatish Balay J[0][1] = ctx->omega_b; 789371c9d4SSatish Balay J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H); 799371c9d4SSatish Balay J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H); 80c4762a1bSJed Brown 819566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 829566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 83c4762a1bSJed Brown 849566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 859566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 86c4762a1bSJed Brown if (A != B) { 879566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 889566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown PetscFunctionReturn(0); 91c4762a1bSJed Brown } 92c4762a1bSJed Brown 93*d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 94*d71ae5a4SJacob Faibussowitsch { 95c4762a1bSJed Brown TS ts; /* ODE integrator */ 96c4762a1bSJed Brown Vec U; /* solution will be stored here */ 97c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 98c4762a1bSJed Brown PetscMPIInt size; 99c4762a1bSJed Brown PetscInt n = 2; 100c4762a1bSJed Brown AppCtx ctx; 101c4762a1bSJed Brown PetscScalar *u; 102c4762a1bSJed Brown PetscReal du[2] = {0.0, 0.0}; 103c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE, flg1, flg2; 104c4762a1bSJed Brown 105c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 106c4762a1bSJed Brown Initialize program 107c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 108327415f7SBarry Smith PetscFunctionBeginUser; 1099566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1109566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1113c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 112c4762a1bSJed Brown 113c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 114c4762a1bSJed Brown Create necessary matrix and vectors 115c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1169566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1179566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1189566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 1199566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1209566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 121c4762a1bSJed Brown 1229566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 123c4762a1bSJed Brown 124c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 125c4762a1bSJed Brown Set runtime options 126c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 127d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 128c4762a1bSJed Brown { 129c4762a1bSJed Brown ctx.omega_b = 1.0; 130c4762a1bSJed Brown ctx.omega_s = 2.0 * PETSC_PI * 60.0; 131c4762a1bSJed Brown ctx.H = 5.0; 1329566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 133c4762a1bSJed Brown ctx.D = 5.0; 1349566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 135c4762a1bSJed Brown ctx.E = 1.1378; 136c4762a1bSJed Brown ctx.V = 1.0; 137c4762a1bSJed Brown ctx.X = 0.545; 138c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 1399566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 140c4762a1bSJed Brown ctx.Pm = 0.9; 1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 142c4762a1bSJed Brown ctx.tf = 1.0; 143c4762a1bSJed Brown ctx.tcl = 1.05; 1449566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 1459566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 1469566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL)); 147c4762a1bSJed Brown if (ensemble) { 148c4762a1bSJed Brown ctx.tf = -1; 149c4762a1bSJed Brown ctx.tcl = -1; 150c4762a1bSJed Brown } 151c4762a1bSJed Brown 1529566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 153c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 154c4762a1bSJed Brown u[1] = 1.0; 1559566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1)); 156c4762a1bSJed Brown n = 2; 1579566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2)); 158c4762a1bSJed Brown u[0] += du[0]; 159c4762a1bSJed Brown u[1] += du[1]; 1609566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 161c4762a1bSJed Brown if (flg1 || flg2) { 162c4762a1bSJed Brown ctx.tf = -1; 163c4762a1bSJed Brown ctx.tcl = -1; 164c4762a1bSJed Brown } 165c4762a1bSJed Brown } 166d0609cedSBarry Smith PetscOptionsEnd(); 167c4762a1bSJed Brown 168c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 169c4762a1bSJed Brown Create timestepping solver context 170c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1719566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1729566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1739566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSTHETA)); 1749566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunction)RHSFunction, &ctx)); 1759566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, (TSRHSJacobian)RHSJacobian, &ctx)); 176c4762a1bSJed Brown 177c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 178c4762a1bSJed Brown Set initial conditions 179c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1809566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 181c4762a1bSJed Brown 182c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 183c4762a1bSJed Brown Set solver options 184c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1859566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 35.0)); 1869566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 1879566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 1889566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 189c4762a1bSJed Brown 190c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 191c4762a1bSJed Brown Solve nonlinear system 192c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 193c4762a1bSJed Brown if (ensemble) { 194c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 1959566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 196c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 197c4762a1bSJed Brown u[1] = ctx.omega_s; 198c4762a1bSJed Brown u[0] += du[0]; 199c4762a1bSJed Brown u[1] += du[1]; 2009566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 2019566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 2029566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 203c4762a1bSJed Brown } 204c4762a1bSJed Brown } else { 2059566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 206c4762a1bSJed Brown } 2079566063dSJacob Faibussowitsch PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD)); 208c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 209c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 210c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2119566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2129566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 2139566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2149566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 215b122ec5aSJacob Faibussowitsch return 0; 216c4762a1bSJed Brown } 217c4762a1bSJed Brown 218c4762a1bSJed Brown /*TEST 219c4762a1bSJed Brown 220c4762a1bSJed Brown build: 221c4762a1bSJed Brown requires: !complex 222c4762a1bSJed Brown 223c4762a1bSJed Brown test: 224c4762a1bSJed Brown 225c4762a1bSJed Brown TEST*/ 226