1c4762a1bSJed Brown static char help[] = "Basic equation for generator stability analysis.\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 Ensemble of initial conditions 11c4762a1bSJed Brown ./ex3 -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 12c4762a1bSJed Brown 13c4762a1bSJed Brown Fault at .1 seconds 14c4762a1bSJed Brown ./ex3 -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 15c4762a1bSJed Brown 16c4762a1bSJed Brown Initial conditions same as when fault is ended 17c4762a1bSJed Brown ./ex3 -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 18c4762a1bSJed Brown 19c4762a1bSJed Brown F*/ 20c4762a1bSJed Brown 21c4762a1bSJed Brown /* 22c4762a1bSJed Brown Include "petscts.h" so that we can use TS solvers. Note that this 23c4762a1bSJed Brown file automatically includes: 24c4762a1bSJed Brown petscsys.h - base PETSc routines petscvec.h - vectors 25c4762a1bSJed Brown petscmat.h - matrices 26c4762a1bSJed Brown petscis.h - index sets petscksp.h - Krylov subspace methods 27c4762a1bSJed Brown petscviewer.h - viewers petscpc.h - preconditioners 28c4762a1bSJed Brown petscksp.h - linear solvers 29c4762a1bSJed Brown */ 30c4762a1bSJed Brown 31c4762a1bSJed Brown #include <petscts.h> 32c4762a1bSJed Brown #include "ex3.h" 33c4762a1bSJed Brown 34d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv) 35d71ae5a4SJacob Faibussowitsch { 36c4762a1bSJed Brown TS ts; /* ODE integrator */ 37c4762a1bSJed Brown Vec U; /* solution will be stored here */ 38c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 39c4762a1bSJed Brown PetscMPIInt size; 40c4762a1bSJed Brown PetscInt n = 2; 41c4762a1bSJed Brown AppCtx ctx; 42c4762a1bSJed Brown PetscScalar *u; 43c4762a1bSJed Brown PetscReal du[2] = {0.0, 0.0}; 44c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE, flg1, flg2; 45c4762a1bSJed Brown PetscInt direction[2]; 46c4762a1bSJed Brown PetscBool terminate[2]; 47c4762a1bSJed Brown 48c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 49c4762a1bSJed Brown Initialize program 50c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 51327415f7SBarry Smith PetscFunctionBeginUser; 529566063dSJacob Faibussowitsch PetscCall(PetscInitialize(&argc, &argv, (char *)0, help)); 539566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 543c633725SBarry Smith PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 55c4762a1bSJed Brown 56c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 57c4762a1bSJed Brown Create necessary matrix and vectors 58c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 599566063dSJacob Faibussowitsch PetscCall(MatCreate(PETSC_COMM_WORLD, &A)); 609566063dSJacob Faibussowitsch PetscCall(MatSetSizes(A, n, n, PETSC_DETERMINE, PETSC_DETERMINE)); 619566063dSJacob Faibussowitsch PetscCall(MatSetType(A, MATDENSE)); 629566063dSJacob Faibussowitsch PetscCall(MatSetFromOptions(A)); 639566063dSJacob Faibussowitsch PetscCall(MatSetUp(A)); 64c4762a1bSJed Brown 659566063dSJacob Faibussowitsch PetscCall(MatCreateVecs(A, &U, NULL)); 66c4762a1bSJed Brown 67c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 68c4762a1bSJed Brown Set runtime options 69c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 70d0609cedSBarry Smith PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Swing equation options", ""); 71c4762a1bSJed Brown { 72c4762a1bSJed Brown ctx.omega_b = 1.0; 73c4762a1bSJed Brown ctx.omega_s = 2.0 * PETSC_PI * 60.0; 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 = 0.9; 859566063dSJacob Faibussowitsch PetscCall(PetscOptionsScalar("-Pm", "", "", ctx.Pm, &ctx.Pm, NULL)); 86c4762a1bSJed Brown ctx.tf = 1.0; 87c4762a1bSJed Brown ctx.tcl = 1.05; 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)); 909566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-ensemble", "Run ensemble of different initial conditions", "", ensemble, &ensemble, NULL)); 91c4762a1bSJed Brown if (ensemble) { 92c4762a1bSJed Brown ctx.tf = -1; 93c4762a1bSJed Brown ctx.tcl = -1; 94c4762a1bSJed Brown } 95c4762a1bSJed Brown 969566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 97c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 98c4762a1bSJed Brown u[1] = 1.0; 999566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-u", "Initial solution", "", u, &n, &flg1)); 100c4762a1bSJed Brown n = 2; 1019566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-du", "Perturbation in initial solution", "", du, &n, &flg2)); 102c4762a1bSJed Brown u[0] += du[0]; 103c4762a1bSJed Brown u[1] += du[1]; 1049566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 105c4762a1bSJed Brown if (flg1 || flg2) { 106c4762a1bSJed Brown ctx.tf = -1; 107c4762a1bSJed Brown ctx.tcl = -1; 108c4762a1bSJed Brown } 109c4762a1bSJed Brown } 110d0609cedSBarry Smith PetscOptionsEnd(); 111c4762a1bSJed Brown 112c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 113c4762a1bSJed Brown Create timestepping solver context 114c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1159566063dSJacob Faibussowitsch PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 1169566063dSJacob Faibussowitsch PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 1179566063dSJacob Faibussowitsch PetscCall(TSSetType(ts, TSTHETA)); 1189566063dSJacob Faibussowitsch PetscCall(TSSetEquationType(ts, TS_EQ_IMPLICIT)); 1199566063dSJacob Faibussowitsch PetscCall(TSARKIMEXSetFullyImplicit(ts, PETSC_TRUE)); 120*8434afd1SBarry Smith PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)IFunction, &ctx)); 121*8434afd1SBarry Smith PetscCall(TSSetIJacobian(ts, A, A, (TSIJacobianFn *)IJacobian, &ctx)); 122c4762a1bSJed Brown 123c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 124c4762a1bSJed Brown Set initial conditions 125c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1269566063dSJacob Faibussowitsch PetscCall(TSSetSolution(ts, U)); 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Set solver options 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1319566063dSJacob Faibussowitsch PetscCall(TSSetMaxTime(ts, 35.0)); 1329566063dSJacob Faibussowitsch PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 1339566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .1)); 1349566063dSJacob Faibussowitsch PetscCall(TSSetFromOptions(ts)); 135c4762a1bSJed Brown 136c4762a1bSJed Brown direction[0] = direction[1] = 1; 137c4762a1bSJed Brown terminate[0] = terminate[1] = PETSC_FALSE; 138c4762a1bSJed Brown 1399566063dSJacob Faibussowitsch PetscCall(TSSetEventHandler(ts, 2, direction, terminate, EventFunction, PostEventFunction, (void *)&ctx)); 140c4762a1bSJed Brown 141c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 142c4762a1bSJed Brown Solve nonlinear system 143c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 144c4762a1bSJed Brown if (ensemble) { 145c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 1469566063dSJacob Faibussowitsch PetscCall(VecGetArray(U, &u)); 147c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm / ctx.Pmax); 148c4762a1bSJed Brown u[1] = ctx.omega_s; 149c4762a1bSJed Brown u[0] += du[0]; 150c4762a1bSJed Brown u[1] += du[1]; 1519566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(U, &u)); 1529566063dSJacob Faibussowitsch PetscCall(TSSetTimeStep(ts, .01)); 1539566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 154c4762a1bSJed Brown } 155c4762a1bSJed Brown } else { 1569566063dSJacob Faibussowitsch PetscCall(TSSolve(ts, U)); 157c4762a1bSJed Brown } 1589566063dSJacob Faibussowitsch PetscCall(VecView(U, PETSC_VIEWER_STDOUT_WORLD)); 159c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 160c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 161c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1629566063dSJacob Faibussowitsch PetscCall(MatDestroy(&A)); 1639566063dSJacob Faibussowitsch PetscCall(VecDestroy(&U)); 1649566063dSJacob Faibussowitsch PetscCall(TSDestroy(&ts)); 1659566063dSJacob Faibussowitsch PetscCall(PetscFinalize()); 166b122ec5aSJacob Faibussowitsch return 0; 167c4762a1bSJed Brown } 168c4762a1bSJed Brown 169c4762a1bSJed Brown /*TEST 170c4762a1bSJed Brown 171c4762a1bSJed Brown build: 172c4762a1bSJed Brown requires: !complex !single 173c4762a1bSJed Brown 174c4762a1bSJed Brown test: 175c4762a1bSJed Brown args: -nox 176c4762a1bSJed Brown 177c4762a1bSJed Brown TEST*/ 178