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 */ 42c4762a1bSJed Brown static PetscErrorCode IFunction(TS ts,PetscReal t,Vec U,Vec Udot,Vec F,AppCtx *ctx) 43c4762a1bSJed Brown { 44c4762a1bSJed Brown PetscErrorCode ierr; 45c4762a1bSJed Brown PetscScalar *f,Pmax; 46c4762a1bSJed Brown const PetscScalar *u,*udot; 47c4762a1bSJed Brown 48c4762a1bSJed Brown PetscFunctionBegin; 49c4762a1bSJed Brown /* The next three lines allow us to access the entries of the vectors directly */ 50c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 51c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 52c4762a1bSJed Brown ierr = VecGetArray(F,&f);CHKERRQ(ierr); 53c4762a1bSJed 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 */ 54c4762a1bSJed Brown else if (t >= ctx->tcl) Pmax = ctx->E/0.745; 55c4762a1bSJed Brown else Pmax = ctx->Pmax; 56c4762a1bSJed Brown f[0] = udot[0] - ctx->omega_s*(u[1] - 1.0); 57c4762a1bSJed Brown f[1] = 2.0*ctx->H*udot[1] + Pmax*PetscSinScalar(u[0]) + ctx->D*(u[1] - 1.0)- ctx->Pm; 58c4762a1bSJed Brown 59c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 60c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 61c4762a1bSJed Brown ierr = VecRestoreArray(F,&f);CHKERRQ(ierr); 62c4762a1bSJed Brown PetscFunctionReturn(0); 63c4762a1bSJed Brown } 64c4762a1bSJed Brown 65c4762a1bSJed Brown /* 66c4762a1bSJed Brown Defines the Jacobian of the ODE passed to the ODE solver. See TSSetIJacobian() for the meaning of a and the Jacobian. 67c4762a1bSJed Brown */ 68c4762a1bSJed Brown static PetscErrorCode IJacobian(TS ts,PetscReal t,Vec U,Vec Udot,PetscReal a,Mat A,Mat B,AppCtx *ctx) 69c4762a1bSJed Brown { 70c4762a1bSJed Brown PetscErrorCode ierr; 71c4762a1bSJed Brown PetscInt rowcol[] = {0,1}; 72c4762a1bSJed Brown PetscScalar J[2][2],Pmax; 73c4762a1bSJed Brown const PetscScalar *u,*udot; 74c4762a1bSJed Brown 75c4762a1bSJed Brown PetscFunctionBegin; 76c4762a1bSJed Brown ierr = VecGetArrayRead(U,&u);CHKERRQ(ierr); 77c4762a1bSJed Brown ierr = VecGetArrayRead(Udot,&udot);CHKERRQ(ierr); 78c4762a1bSJed 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 */ 79c4762a1bSJed Brown else if (t >= ctx->tcl) Pmax = ctx->E/0.745; 80c4762a1bSJed Brown else Pmax = ctx->Pmax; 81c4762a1bSJed Brown 82c4762a1bSJed Brown J[0][0] = a; J[0][1] = -ctx->omega_s; 83c4762a1bSJed Brown J[1][1] = 2.0*ctx->H*a + ctx->D; J[1][0] = Pmax*PetscCosScalar(u[0]); 84c4762a1bSJed Brown 85c4762a1bSJed Brown ierr = MatSetValues(B,2,rowcol,2,rowcol,&J[0][0],INSERT_VALUES);CHKERRQ(ierr); 86c4762a1bSJed Brown ierr = VecRestoreArrayRead(U,&u);CHKERRQ(ierr); 87c4762a1bSJed Brown ierr = VecRestoreArrayRead(Udot,&udot);CHKERRQ(ierr); 88c4762a1bSJed Brown 89c4762a1bSJed Brown ierr = MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 90c4762a1bSJed Brown ierr = MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 91c4762a1bSJed Brown if (A != B) { 92c4762a1bSJed Brown ierr = MatAssemblyBegin(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 93c4762a1bSJed Brown ierr = MatAssemblyEnd(B,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 94c4762a1bSJed Brown } 95c4762a1bSJed Brown PetscFunctionReturn(0); 96c4762a1bSJed Brown } 97c4762a1bSJed Brown 98c4762a1bSJed Brown PetscErrorCode PostStep(TS ts) 99c4762a1bSJed Brown { 100c4762a1bSJed Brown PetscErrorCode ierr; 101c4762a1bSJed Brown Vec X; 102c4762a1bSJed Brown PetscReal t; 103c4762a1bSJed Brown 104c4762a1bSJed Brown PetscFunctionBegin; 105c4762a1bSJed Brown ierr = TSGetTime(ts,&t);CHKERRQ(ierr); 106c4762a1bSJed Brown if (t >= .2) { 107c4762a1bSJed Brown ierr = TSGetSolution(ts,&X);CHKERRQ(ierr); 108c4762a1bSJed Brown ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 109c4762a1bSJed Brown exit(0); 110c4762a1bSJed Brown /* results in initial conditions after fault of -u 0.496792,1.00932 */ 111c4762a1bSJed Brown } 112c4762a1bSJed Brown PetscFunctionReturn(0); 113c4762a1bSJed Brown } 114c4762a1bSJed Brown 115c4762a1bSJed Brown int main(int argc,char **argv) 116c4762a1bSJed Brown { 117c4762a1bSJed Brown TS ts; /* ODE integrator */ 118c4762a1bSJed Brown Vec U; /* solution will be stored here */ 119c4762a1bSJed Brown Mat A; /* Jacobian matrix */ 120c4762a1bSJed Brown PetscErrorCode ierr; 121c4762a1bSJed Brown PetscMPIInt size; 122c4762a1bSJed Brown PetscInt n = 2; 123c4762a1bSJed Brown AppCtx ctx; 124c4762a1bSJed Brown PetscScalar *u; 125c4762a1bSJed Brown PetscReal du[2] = {0.0,0.0}; 126c4762a1bSJed Brown PetscBool ensemble = PETSC_FALSE,flg1,flg2; 127c4762a1bSJed Brown 128c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 129c4762a1bSJed Brown Initialize program 130c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 131c4762a1bSJed Brown ierr = PetscInitialize(&argc,&argv,(char*)0,help);if (ierr) return ierr; 132ffc4695bSBarry Smith ierr = MPI_Comm_size(PETSC_COMM_WORLD,&size);CHKERRMPI(ierr); 133*3c633725SBarry Smith PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 134c4762a1bSJed Brown 135c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 136c4762a1bSJed Brown Create necessary matrix and vectors 137c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 138c4762a1bSJed Brown ierr = MatCreate(PETSC_COMM_WORLD,&A);CHKERRQ(ierr); 139c4762a1bSJed Brown ierr = MatSetSizes(A,n,n,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 140c4762a1bSJed Brown ierr = MatSetType(A,MATDENSE);CHKERRQ(ierr); 141c4762a1bSJed Brown ierr = MatSetFromOptions(A);CHKERRQ(ierr); 142c4762a1bSJed Brown ierr = MatSetUp(A);CHKERRQ(ierr); 143c4762a1bSJed Brown 144c4762a1bSJed Brown ierr = MatCreateVecs(A,&U,NULL);CHKERRQ(ierr); 145c4762a1bSJed Brown 146c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 147c4762a1bSJed Brown Set runtime options 148c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 149c4762a1bSJed Brown ierr = PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Swing equation options","");CHKERRQ(ierr); 150c4762a1bSJed Brown { 151c4762a1bSJed Brown ctx.omega_s = 2.0*PETSC_PI*60.0; 152c4762a1bSJed Brown ctx.H = 5.0; 153c4762a1bSJed Brown ierr = PetscOptionsScalar("-Inertia","","",ctx.H,&ctx.H,NULL);CHKERRQ(ierr); 154c4762a1bSJed Brown ctx.D = 5.0; 155c4762a1bSJed Brown ierr = PetscOptionsScalar("-D","","",ctx.D,&ctx.D,NULL);CHKERRQ(ierr); 156c4762a1bSJed Brown ctx.E = 1.1378; 157c4762a1bSJed Brown ctx.V = 1.0; 158c4762a1bSJed Brown ctx.X = 0.545; 159c4762a1bSJed Brown ctx.Pmax = ctx.E*ctx.V/ctx.X; 160c4762a1bSJed Brown ierr = PetscOptionsScalar("-Pmax","","",ctx.Pmax,&ctx.Pmax,NULL);CHKERRQ(ierr); 161c4762a1bSJed Brown ctx.Pm = 0.9; 162c4762a1bSJed Brown ierr = PetscOptionsScalar("-Pm","","",ctx.Pm,&ctx.Pm,NULL);CHKERRQ(ierr); 163c4762a1bSJed Brown ctx.tf = 1.0; 164c4762a1bSJed Brown ctx.tcl = 1.05; 165c4762a1bSJed Brown ierr = PetscOptionsReal("-tf","Time to start fault","",ctx.tf,&ctx.tf,NULL);CHKERRQ(ierr); 166c4762a1bSJed Brown ierr = PetscOptionsReal("-tcl","Time to end fault","",ctx.tcl,&ctx.tcl,NULL);CHKERRQ(ierr); 167c4762a1bSJed Brown ierr = PetscOptionsBool("-ensemble","Run ensemble of different initial conditions","",ensemble,&ensemble,NULL);CHKERRQ(ierr); 168c4762a1bSJed Brown if (ensemble) { 169c4762a1bSJed Brown ctx.tf = -1; 170c4762a1bSJed Brown ctx.tcl = -1; 171c4762a1bSJed Brown } 172c4762a1bSJed Brown 173c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 174c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 175c4762a1bSJed Brown u[1] = 1.0; 176c4762a1bSJed Brown ierr = PetscOptionsRealArray("-u","Initial solution","",u,&n,&flg1);CHKERRQ(ierr); 177c4762a1bSJed Brown n = 2; 178c4762a1bSJed Brown ierr = PetscOptionsRealArray("-du","Perturbation in initial solution","",du,&n,&flg2);CHKERRQ(ierr); 179c4762a1bSJed Brown u[0] += du[0]; 180c4762a1bSJed Brown u[1] += du[1]; 181c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 182c4762a1bSJed Brown if (flg1 || flg2) { 183c4762a1bSJed Brown ctx.tf = -1; 184c4762a1bSJed Brown ctx.tcl = -1; 185c4762a1bSJed Brown } 186c4762a1bSJed Brown } 187c4762a1bSJed Brown ierr = PetscOptionsEnd();CHKERRQ(ierr); 188c4762a1bSJed Brown 189c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 190c4762a1bSJed Brown Create timestepping solver context 191c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 192c4762a1bSJed Brown ierr = TSCreate(PETSC_COMM_WORLD,&ts);CHKERRQ(ierr); 193c4762a1bSJed Brown ierr = TSSetProblemType(ts,TS_NONLINEAR);CHKERRQ(ierr); 194c4762a1bSJed Brown ierr = TSSetType(ts,TSROSW);CHKERRQ(ierr); 195c4762a1bSJed Brown ierr = TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&ctx);CHKERRQ(ierr); 196c4762a1bSJed Brown ierr = TSSetIJacobian(ts,A,A,(TSIJacobian)IJacobian,&ctx);CHKERRQ(ierr); 197c4762a1bSJed Brown 198c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 199c4762a1bSJed Brown Set initial conditions 200c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 201c4762a1bSJed Brown ierr = TSSetSolution(ts,U);CHKERRQ(ierr); 202c4762a1bSJed Brown 203c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 204c4762a1bSJed Brown Set solver options 205c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 206c4762a1bSJed Brown ierr = TSSetMaxTime(ts,35.0);CHKERRQ(ierr); 207c4762a1bSJed Brown ierr = TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP);CHKERRQ(ierr); 208c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr); 209c4762a1bSJed Brown ierr = TSSetFromOptions(ts);CHKERRQ(ierr); 210c4762a1bSJed Brown /* ierr = TSSetPostStep(ts,PostStep);CHKERRQ(ierr); */ 211c4762a1bSJed Brown 212c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 213c4762a1bSJed Brown Solve nonlinear system 214c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 215c4762a1bSJed Brown if (ensemble) { 216c4762a1bSJed Brown for (du[1] = -2.5; du[1] <= .01; du[1] += .1) { 217c4762a1bSJed Brown ierr = VecGetArray(U,&u);CHKERRQ(ierr); 218c4762a1bSJed Brown u[0] = PetscAsinScalar(ctx.Pm/ctx.Pmax); 219c4762a1bSJed Brown u[1] = ctx.omega_s; 220c4762a1bSJed Brown u[0] += du[0]; 221c4762a1bSJed Brown u[1] += du[1]; 222c4762a1bSJed Brown ierr = VecRestoreArray(U,&u);CHKERRQ(ierr); 223c4762a1bSJed Brown ierr = TSSetTimeStep(ts,.01);CHKERRQ(ierr); 224c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 225c4762a1bSJed Brown } 226c4762a1bSJed Brown } else { 227c4762a1bSJed Brown ierr = TSSolve(ts,U);CHKERRQ(ierr); 228c4762a1bSJed Brown } 229c4762a1bSJed Brown /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 230c4762a1bSJed Brown Free work space. All PETSc objects should be destroyed when they are no longer needed. 231c4762a1bSJed Brown - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 232c4762a1bSJed Brown ierr = MatDestroy(&A);CHKERRQ(ierr); 233c4762a1bSJed Brown ierr = VecDestroy(&U);CHKERRQ(ierr); 234c4762a1bSJed Brown ierr = TSDestroy(&ts);CHKERRQ(ierr); 235c4762a1bSJed Brown ierr = PetscFinalize(); 236c4762a1bSJed Brown return ierr; 237c4762a1bSJed Brown } 238c4762a1bSJed Brown 239c4762a1bSJed Brown /*TEST 240c4762a1bSJed Brown 241c4762a1bSJed Brown build: 242c4762a1bSJed Brown requires: !complex 243c4762a1bSJed Brown 244c4762a1bSJed Brown test: 245c4762a1bSJed Brown args: -nox -ts_dt 10 246c4762a1bSJed Brown 247c4762a1bSJed Brown TEST*/ 248