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