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 */ 42c4762a1bSJed Brown static PetscErrorCode RHSFunction(TS ts,PetscReal t,Vec U,Vec F,AppCtx *ctx) 43c4762a1bSJed Brown { 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 */ 495f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArrayRead(U,&u)); 505f80ce2aSJacob Faibussowitsch CHKERRQ(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 575f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 585f80ce2aSJacob Faibussowitsch CHKERRQ(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 */ 65c4762a1bSJed Brown static PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec U,Mat A,Mat B,AppCtx *ctx) 66c4762a1bSJed Brown { 67c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 68c4762a1bSJed Brown PetscScalar J[2][2],Pmax; 69c4762a1bSJed Brown const PetscScalar *u; 70c4762a1bSJed Brown 71c4762a1bSJed Brown PetscFunctionBegin; 725f80ce2aSJacob Faibussowitsch CHKERRQ(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 76c4762a1bSJed Brown J[0][0] = 0; J[0][1] = ctx->omega_b; 77c4762a1bSJed Brown J[1][1] = -ctx->D*ctx->omega_s/(2.0*ctx->H); J[1][0] = -Pmax*PetscCosScalar(u[0])*ctx->omega_s/(2.0*ctx->H); 78c4762a1bSJed Brown 795f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES)); 805f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArrayRead(U,&u)); 81c4762a1bSJed Brown 825f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 835f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 84c4762a1bSJed Brown if (A != B) { 855f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY)); 865f80ce2aSJacob Faibussowitsch CHKERRQ(MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY)); 87c4762a1bSJed Brown } 88c4762a1bSJed Brown PetscFunctionReturn(0); 89c4762a1bSJed Brown } 90c4762a1bSJed Brown 91c4762a1bSJed Brown int main(int argc,char **argv) 92c4762a1bSJed Brown { 93c4762a1bSJed Brown TS ts; /* ODE integrator */ 94c4762a1bSJed Brown Vec U; /* solution will be stored here */ 95c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 96c4762a1bSJed Brown PetscErrorCode ierr; 97c4762a1bSJed Brown PetscMPIInt size; 98c4762a1bSJed Brown PetscInt n = 2; 99c4762a1bSJed Brown AppCtx ctx; 100c4762a1bSJed Brown PetscScalar *u; 101c4762a1bSJed Brown PetscReal du[2] = {0.0,0.0}; 102c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE,flg1,flg2; 103c4762a1bSJed Brown 104c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 105c4762a1bSJed Brown Initialize program 106c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 107*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscInitialize(&argc,&argv,(char*)0,help)); 1085f80ce2aSJacob Faibussowitsch CHKERRMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1093c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 110c4762a1bSJed Brown 111c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 112c4762a1bSJed Brown Create necessary matrix and vectors 113c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1145f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreate(PETSC_COMM_WORLD,&A)); 1155f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE)); 1165f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetType(A,MATDENSE)); 1175f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetFromOptions(A)); 1185f80ce2aSJacob Faibussowitsch CHKERRQ(MatSetUp(A)); 119c4762a1bSJed Brown 1205f80ce2aSJacob Faibussowitsch CHKERRQ(MatCreateVecs(A,&U,NULL)); 121c4762a1bSJed Brown 122c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 123c4762a1bSJed Brown Set runtime options 124c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 125c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr); 126c4762a1bSJed Brown { 127c4762a1bSJed Brown ctx.omega_b = 1.0; 128c4762a1bSJed Brown ctx.omega_s = 2.0*PETSC_PI*60.0; 129c4762a1bSJed Brown ctx.H = 5.0; 1305f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL)); 131c4762a1bSJed Brown ctx.D = 5.0; 1325f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL)); 133c4762a1bSJed Brown ctx.E = 1.1378; 134c4762a1bSJed Brown ctx.V = 1.0; 135c4762a1bSJed Brown ctx.X = 0.545; 136c4762a1bSJed Brown ctx.Pmax = ctx.E*ctx.V/ctx.X; 1375f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL)); 138c4762a1bSJed Brown ctx.Pm = 0.9; 1395f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL)); 140c4762a1bSJed Brown ctx.tf = 1.0; 141c4762a1bSJed Brown ctx.tcl = 1.05; 1425f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL)); 1435f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL)); 1445f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL)); 145c4762a1bSJed Brown if (ensemble) { 146c4762a1bSJed Brown ctx.tf = -1; 147c4762a1bSJed Brown ctx.tcl = -1; 148c4762a1bSJed Brown } 149c4762a1bSJed Brown 1505f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(U,&u)); 151c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 152c4762a1bSJed Brown u[1] = 1.0; 1535f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1)); 154c4762a1bSJed Brown n = 2; 1555f80ce2aSJacob Faibussowitsch CHKERRQ(PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2)); 156c4762a1bSJed Brown u[0] += du[0]; 157c4762a1bSJed Brown u[1] += du[1]; 1585f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(U,&u)); 159c4762a1bSJed Brown if (flg1 || flg2) { 160c4762a1bSJed Brown ctx.tf = -1; 161c4762a1bSJed Brown ctx.tcl = -1; 162c4762a1bSJed Brown } 163c4762a1bSJed Brown } 164c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 165c4762a1bSJed Brown 166c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 167c4762a1bSJed Brown Create timestepping solver context 168c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1695f80ce2aSJacob Faibussowitsch CHKERRQ(TSCreate(PETSC_COMM_WORLD,&ts)); 1705f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetProblemType(ts,TS_NONLINEAR)); 1715f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetType(ts,TSTHETA)); 1725f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSFunction(ts,NULL,(TSRHSFunction)RHSFunction,&ctx)); 1735f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetRHSJacobian(ts,A,A,(TSRHSJacobian)RHSJacobian,&ctx)); 174c4762a1bSJed Brown 175c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 176c4762a1bSJed Brown Set initial conditions 177c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1785f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetSolution(ts,U)); 179c4762a1bSJed Brown 180c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 181c4762a1bSJed Brown Set solver options 182c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1835f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetMaxTime(ts,35.0)); 1845f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 1855f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,.01)); 1865f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetFromOptions(ts)); 187c4762a1bSJed Brown 188c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 189c4762a1bSJed Brown Solve nonlinear system 190c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 191c4762a1bSJed Brown if (ensemble) { 192c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 1935f80ce2aSJacob Faibussowitsch CHKERRQ(VecGetArray(U,&u)); 194c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 195c4762a1bSJed Brown u[1] = ctx.omega_s; 196c4762a1bSJed Brown u[0] += du[0]; 197c4762a1bSJed Brown u[1] += du[1]; 1985f80ce2aSJacob Faibussowitsch CHKERRQ(VecRestoreArray(U,&u)); 1995f80ce2aSJacob Faibussowitsch CHKERRQ(TSSetTimeStep(ts,.01)); 2005f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,U)); 201c4762a1bSJed Brown } 202c4762a1bSJed Brown } else { 2035f80ce2aSJacob Faibussowitsch CHKERRQ(TSSolve(ts,U)); 204c4762a1bSJed Brown } 2055f80ce2aSJacob Faibussowitsch CHKERRQ(VecView(U,PETSC_VIEWER_STDOUT_WORLD)); 206c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 207c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 208c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 2095f80ce2aSJacob Faibussowitsch CHKERRQ(MatDestroy(&A)); 2105f80ce2aSJacob Faibussowitsch CHKERRQ(VecDestroy(&U)); 2115f80ce2aSJacob Faibussowitsch CHKERRQ(TSDestroy(&ts)); 212*b122ec5aSJacob Faibussowitsch CHKERRQ(PetscFinalize()); 213*b122ec5aSJacob Faibussowitsch return 0; 214c4762a1bSJed Brown } 215c4762a1bSJed Brown 216c4762a1bSJed Brown /*TEST 217c4762a1bSJed Brown 218c4762a1bSJed Brown build: 219c4762a1bSJed Brown requires: !complex 220c4762a1bSJed Brown 221c4762a1bSJed Brown test: 222c4762a1bSJed Brown 223c4762a1bSJed Brown TEST*/ 224