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*9371c9d4SSatish Balay static PetscErrorCode RHSFunction(TS ts, PetscReal t, Vec U, Vec F, AppCtx *ctx) { 43c4762a1bSJed Brown const PetscScalar *u; 44c4762a1bSJed Brown PetscScalar *f, Pmax; 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(VecGetArray(F, &f)); 50c4762a1bSJed 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 */ 51c4762a1bSJed Brown else Pmax = ctx->Pmax; 52c4762a1bSJed Brown 53c4762a1bSJed Brown f[0] = ctx->omega_b * (u[1] - ctx->omega_s); 54c4762a1bSJed Brown f[1] = (-Pmax * PetscSinScalar(u[0]) - ctx->D * (u[1] - ctx->omega_s) + ctx->Pm) * ctx->omega_s / (2.0 * ctx->H); 55c4762a1bSJed Brown 569566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(F, &f)); 58c4762a1bSJed Brown PetscFunctionReturn(0); 59c4762a1bSJed Brown } 60c4762a1bSJed Brown 61c4762a1bSJed Brown /* 62c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 63c4762a1bSJed Brown */ 64*9371c9d4SSatish Balay static PetscErrorCode RHSJacobian(TS ts, PetscReal t, Vec U, Mat A, Mat B, AppCtx *ctx) { 65c4762a1bSJed Brown PetscInt rowcol[] = {0, 1}; 66c4762a1bSJed Brown PetscScalar J[2][2], Pmax; 67c4762a1bSJed Brown const PetscScalar *u; 68c4762a1bSJed Brown 69c4762a1bSJed Brown PetscFunctionBegin; 709566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(U, &u)); 71c4762a1bSJed 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 */ 72c4762a1bSJed Brown else Pmax = ctx->Pmax; 73c4762a1bSJed Brown 74*9371c9d4SSatish Balay J[0][0] = 0; 75*9371c9d4SSatish Balay J[0][1] = ctx->omega_b; 76*9371c9d4SSatish Balay J[1][1] = -ctx->D * ctx->omega_s / (2.0 * ctx->H); 77*9371c9d4SSatish Balay J[1][0] = -Pmax * PetscCosScalar(u[0]) * ctx->omega_s / (2.0 * ctx->H); 78c4762a1bSJed Brown 799566063dSJacob Faibussowitsch PetscCall(MatSetValues(B, 2, rowcol, 2, rowcol, &J[0][0], INSERT_VALUES)); 809566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(U, &u)); 81c4762a1bSJed Brown 829566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 839566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 84c4762a1bSJed Brown if (A != B) { 859566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(B, MAT_FINAL_ASSEMBLY)); 869566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(B, MAT_FINAL_ASSEMBLY)); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown PetscFunctionReturn(0); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91*9371c9d4SSatish Balay int main(int argc, char **argv) { 92c4762a1bSJed Brown TS ts; /* ODE integrator */ 93c4762a1bSJed Brown Vec U; /* solution will be stored here */ 94c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 95c4762a1bSJed Brown PetscMPIInt size; 96c4762a1bSJed Brown PetscInt n = 2; 97c4762a1bSJed Brown AppCtx ctx; 98c4762a1bSJed Brown PetscScalar *u; 99c4762a1bSJed Brown PetscReal du[2] = {0.0, 0.0}; 100c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE, flg1, flg2; 101c4762a1bSJed Brown 102c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 103c4762a1bSJed Brown Initialize program 104c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 105327415f7SBarry Smith PetscFunctionBeginUser; 1069566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 1079566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 1083c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 109c4762a1bSJed Brown 110c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 111c4762a1bSJed Brown Create necessary matrix and vectors 112c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1139566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 1149566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 1159566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 1169566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 1179566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 118c4762a1bSJed Brown 1199566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 120c4762a1bSJed Brown 121c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 122c4762a1bSJed Brown Set runtime options 123c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 124d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 125c4762a1bSJed Brown { 126c4762a1bSJed Brown ctx.omega_b = 1.0; 127c4762a1bSJed Brown ctx.omega_s = 2.0 * PETSC_PI * 60.0; 128c4762a1bSJed Brown ctx.H = 5.0; 1299566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Inertia", "", "", ctx.H, &ctx.H, NULL)); 130c4762a1bSJed Brown ctx.D = 5.0; 1319566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-D", "", "", ctx.D, &ctx.D, NULL)); 132c4762a1bSJed Brown ctx.E = 1.1378; 133c4762a1bSJed Brown ctx.V = 1.0; 134c4762a1bSJed Brown ctx.X = 0.545; 135c4762a1bSJed Brown ctx.Pmax = ctx.E * ctx.V / ctx.X; 1369566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pmax", "", "", ctx.Pmax, &ctx.Pmax, NULL)); 137c4762a1bSJed Brown ctx.Pm = 0.9; 1389566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 139c4762a1bSJed Brown ctx.tf = 1.0; 140c4762a1bSJed Brown ctx.tcl = 1.05; 1419566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tf", "Time to start fault", "", ctx.tf, &ctx.tf, NULL)); 1429566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-tcl", "Time to end fault", "", ctx.tcl, &ctx.tcl, NULL)); 1439566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL)); 144c4762a1bSJed Brown if (ensemble) { 145c4762a1bSJed Brown ctx.tf = -1; 146c4762a1bSJed Brown ctx.tcl = -1; 147c4762a1bSJed Brown } 148c4762a1bSJed Brown 1499566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 150c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 151c4762a1bSJed Brown u[1] = 1.0; 1529566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1)); 153c4762a1bSJed Brown n = 2; 1549566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2)); 155c4762a1bSJed Brown u[0] += du[0]; 156c4762a1bSJed Brown u[1] += du[1]; 1579566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 158c4762a1bSJed Brown if (flg1 || flg2) { 159c4762a1bSJed Brown ctx.tf = -1; 160c4762a1bSJed Brown ctx.tcl = -1; 161c4762a1bSJed Brown } 162c4762a1bSJed Brown } 163d0609cedSBarry Smith PetscOptionsEnd(); 164c4762a1bSJed Brown 165c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 166c4762a1bSJed Brown Create timestepping solver context 167c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1689566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1699566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1709566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSTHETA)); 1719566063dSJacob Faibussowitsch PetscCall(TSSetRHSFunction(ts, NULL, (TSRHSFunction)RHSFunction, &ctx)); 1729566063dSJacob Faibussowitsch PetscCall(TSSetRHSJacobian(ts, A, A, (TSRHSJacobian)RHSJacobian, &ctx)); 173c4762a1bSJed Brown 174c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 175c4762a1bSJed Brown Set initial conditions 176c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1779566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 178c4762a1bSJed Brown 179c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 180c4762a1bSJed Brown Set solver options 181c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1829566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 35.0)); 1839566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 1849566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 1859566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 186c4762a1bSJed Brown 187c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 188c4762a1bSJed Brown Solve nonlinear system 189c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 190c4762a1bSJed Brown if (ensemble) { 191c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 1929566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 193c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 194c4762a1bSJed Brown u[1] = ctx.omega_s; 195c4762a1bSJed Brown u[0] += du[0]; 196c4762a1bSJed Brown u[1] += du[1]; 1979566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 1989566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 1999566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 200c4762a1bSJed Brown } 201c4762a1bSJed Brown } else { 2029566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 203c4762a1bSJed Brown } 2049566063dSJacob Faibussowitsch PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD)); 205c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 206c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 207c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2089566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 2099566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 2109566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 2119566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 212b122ec5aSJacob Faibussowitsch return 0; 213c4762a1bSJed Brown } 214c4762a1bSJed Brown 215c4762a1bSJed Brown /*TEST 216c4762a1bSJed Brown 217c4762a1bSJed Brown build: 218c4762a1bSJed Brown requires: !complex 219c4762a1bSJed Brown 220c4762a1bSJed Brown test: 221c4762a1bSJed Brown 222c4762a1bSJed Brown TEST*/ 223