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