1 2 static char help[] = "Power grid stability analysis of WECC 9 bus system.\n\ 3 This example is based on the 9-bus (node) example given in the book Power\n\ 4 Systems Dynamics and Stability (Chapter 7) by P. Sauer and M. A. Pai.\n\ 5 The power grid in this example consists of 9 buses (nodes), 3 generators,\n\ 6 3 loads, and 9 transmission lines. The network equations are written\n\ 7 in current balance form using rectangular coordinates.\n\n"; 8 9 /* 10 The equations for the stability analysis are described by the DAE 11 12 \dot{x} = f(x,y,t) 13 0 = g(x,y,t) 14 15 where the generators are described by differential equations, while the algebraic 16 constraints define the network equations. 17 18 The generators are modeled with a 4th order differential equation describing the electrical 19 and mechanical dynamics. Each generator also has an exciter system modeled by 3rd order 20 diff. eqns. describing the exciter, voltage regulator, and the feedback stabilizer 21 mechanism. 22 23 The network equations are described by nodal current balance equations. 24 I(x,y) - Y*V = 0 25 26 where: 27 I(x,y) is the current injected from generators and loads. 28 Y is the admittance matrix, and 29 V is the voltage vector 30 */ 31 32 /* 33 Include "petscts.h" so that we can use TS solvers. Note that this 34 file automatically includes: 35 petscsys.h - base PETSc routines petscvec.h - vectors 36 petscmat.h - matrices 37 petscis.h - index sets petscksp.h - Krylov subspace methods 38 petscviewer.h - viewers petscpc.h - preconditioners 39 petscksp.h - linear solvers 40 */ 41 42 #include <petscts.h> 43 #include <petscdm.h> 44 #include <petscdmda.h> 45 #include <petscdmcomposite.h> 46 47 #define freq 60 48 #define w_s (2*PETSC_PI*freq) 49 50 /* Sizes and indices */ 51 const PetscInt nbus = 9; /* Number of network buses */ 52 const PetscInt ngen = 3; /* Number of generators */ 53 const PetscInt nload = 3; /* Number of loads */ 54 const PetscInt gbus[3] = {0,1,2}; /* Buses at which generators are incident */ 55 const PetscInt lbus[3] = {4,5,7}; /* Buses at which loads are incident */ 56 57 /* Generator real and reactive powers (found via loadflow) */ 58 const PetscScalar PG[3] = {0.716786142395021,1.630000000000000,0.850000000000000}; 59 const PetscScalar QG[3] = {0.270702180178785,0.066120127797275,-0.108402221791588}; 60 /* Generator constants */ 61 const PetscScalar H[3] = {23.64,6.4,3.01}; /* Inertia constant */ 62 const PetscScalar Rs[3] = {0.0,0.0,0.0}; /* Stator Resistance */ 63 const PetscScalar Xd[3] = {0.146,0.8958,1.3125}; /* d-axis reactance */ 64 const PetscScalar Xdp[3] = {0.0608,0.1198,0.1813}; /* d-axis transient reactance */ 65 const PetscScalar Xq[3] = {0.4360,0.8645,1.2578}; /* q-axis reactance Xq(1) set to 0.4360, value given in text 0.0969 */ 66 const PetscScalar Xqp[3] = {0.0969,0.1969,0.25}; /* q-axis transient reactance */ 67 const PetscScalar Td0p[3] = {8.96,6.0,5.89}; /* d-axis open circuit time constant */ 68 const PetscScalar Tq0p[3] = {0.31,0.535,0.6}; /* q-axis open circuit time constant */ 69 PetscScalar M[3]; /* M = 2*H/w_s */ 70 PetscScalar D[3]; /* D = 0.1*M */ 71 72 PetscScalar TM[3]; /* Mechanical Torque */ 73 /* Exciter system constants */ 74 const PetscScalar KA[3] = {20.0,20.0,20.0}; /* Voltage regulartor gain constant */ 75 const PetscScalar TA[3] = {0.2,0.2,0.2}; /* Voltage regulator time constant */ 76 const PetscScalar KE[3] = {1.0,1.0,1.0}; /* Exciter gain constant */ 77 const PetscScalar TE[3] = {0.314,0.314,0.314}; /* Exciter time constant */ 78 const PetscScalar KF[3] = {0.063,0.063,0.063}; /* Feedback stabilizer gain constant */ 79 const PetscScalar TF[3] = {0.35,0.35,0.35}; /* Feedback stabilizer time constant */ 80 const PetscScalar k1[3] = {0.0039,0.0039,0.0039}; 81 const PetscScalar k2[3] = {1.555,1.555,1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */ 82 const PetscScalar VRMIN[3] = {-4.0,-4.0,-4.0}; 83 const PetscScalar VRMAX[3] = {7.0,7.0,7.0}; 84 PetscInt VRatmin[3]; 85 PetscInt VRatmax[3]; 86 87 PetscScalar Vref[3]; 88 /* Load constants 89 We use a composite load model that describes the load and reactive powers at each time instant as follows 90 P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i 91 Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i 92 where 93 ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads 94 ld_alphap,ld_alphap - Percentage contribution (weights) or loads 95 P_D0 - Real power load 96 Q_D0 - Reactive power load 97 V_m(t) - Voltage magnitude at time t 98 V_m0 - Voltage magnitude at t = 0 99 ld_betap, ld_betaq - exponents describing the load model for real and reactive part 100 101 Note: All loads have the same characteristic currently. 102 */ 103 const PetscScalar PD0[3] = {1.25,0.9,1.0}; 104 const PetscScalar QD0[3] = {0.5,0.3,0.35}; 105 const PetscInt ld_nsegsp[3] = {3,3,3}; 106 const PetscScalar ld_alphap[3] = {1.0,0.0,0.0}; 107 const PetscScalar ld_betap[3] = {2.0,1.0,0.0}; 108 const PetscInt ld_nsegsq[3] = {3,3,3}; 109 const PetscScalar ld_alphaq[3] = {1.0,0.0,0.0}; 110 const PetscScalar ld_betaq[3] = {2.0,1.0,0.0}; 111 112 typedef struct { 113 DM dmgen, dmnet; /* DMs to manage generator and network subsystem */ 114 DM dmpgrid; /* Composite DM to manage the entire power grid */ 115 Mat Ybus; /* Network admittance matrix */ 116 Vec V0; /* Initial voltage vector (Power flow solution) */ 117 PetscReal tfaulton,tfaultoff; /* Fault on and off times */ 118 PetscInt faultbus; /* Fault bus */ 119 PetscScalar Rfault; 120 PetscReal t0,tmax; 121 PetscInt neqs_gen,neqs_net,neqs_pgrid; 122 Mat Sol; /* Matrix to save solution at each time step */ 123 PetscInt stepnum; 124 PetscReal t; 125 SNES snes_alg; 126 IS is_diff; /* indices for differential equations */ 127 IS is_alg; /* indices for algebraic equations */ 128 PetscBool setisdiff; /* TS computes truncation error based only on the differential variables */ 129 PetscBool semiexplicit; /* If the flag is set then a semi-explicit method is used using TSRK */ 130 } Userctx; 131 132 /* 133 The first two events are for fault on and off, respectively. The following events are 134 to check the min/max limits on the state variable VR. A non windup limiter is used for 135 the VR limits. 136 */ 137 PetscErrorCode EventFunction(TS ts,PetscReal t,Vec X,PetscScalar *fvalue,void *ctx) 138 { 139 Userctx *user=(Userctx*)ctx; 140 Vec Xgen,Xnet; 141 PetscInt i,idx=0; 142 const PetscScalar *xgen,*xnet; 143 PetscScalar Efd,RF,VR,Vr,Vi,Vm; 144 145 PetscFunctionBegin; 146 147 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 148 PetscCall(DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet)); 149 150 PetscCall(VecGetArrayRead(Xgen,&xgen)); 151 PetscCall(VecGetArrayRead(Xnet,&xnet)); 152 153 /* Event for fault-on time */ 154 fvalue[0] = t - user->tfaulton; 155 /* Event for fault-off time */ 156 fvalue[1] = t - user->tfaultoff; 157 158 for (i=0; i < ngen; i++) { 159 Efd = xgen[idx+6]; 160 RF = xgen[idx+7]; 161 VR = xgen[idx+8]; 162 163 Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ 164 Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ 165 Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); 166 167 if (!VRatmax[i]) { 168 fvalue[2+2*i] = VRMAX[i] - VR; 169 } else { 170 fvalue[2+2*i] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; 171 } 172 if (!VRatmin[i]) { 173 fvalue[2+2*i+1] = VRMIN[i] - VR; 174 } else { 175 fvalue[2+2*i+1] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; 176 } 177 idx = idx+9; 178 } 179 PetscCall(VecRestoreArrayRead(Xgen,&xgen)); 180 PetscCall(VecRestoreArrayRead(Xnet,&xnet)); 181 182 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 183 184 PetscFunctionReturn(0); 185 } 186 187 PetscErrorCode PostEventFunction(TS ts,PetscInt nevents,PetscInt event_list[],PetscReal t,Vec X,PetscBool forwardsolve,void* ctx) 188 { 189 Userctx *user=(Userctx*)ctx; 190 Vec Xgen,Xnet; 191 PetscScalar *xgen,*xnet; 192 PetscInt row_loc,col_loc; 193 PetscScalar val; 194 PetscInt i,idx=0,event_num; 195 PetscScalar fvalue; 196 PetscScalar Efd, RF, VR; 197 PetscScalar Vr,Vi,Vm; 198 199 PetscFunctionBegin; 200 201 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 202 PetscCall(DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet)); 203 204 PetscCall(VecGetArray(Xgen,&xgen)); 205 PetscCall(VecGetArray(Xnet,&xnet)); 206 207 for (i=0; i < nevents; i++) { 208 if (event_list[i] == 0) { 209 /* Apply disturbance - resistive fault at user->faultbus */ 210 /* This is done by adding shunt conductance to the diagonal location 211 in the Ybus matrix */ 212 row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1; /* Location for G */ 213 val = 1/user->Rfault; 214 PetscCall(MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES)); 215 row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus; /* Location for G */ 216 val = 1/user->Rfault; 217 PetscCall(MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES)); 218 219 PetscCall(MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY)); 220 PetscCall(MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY)); 221 222 /* Solve the algebraic equations */ 223 PetscCall(SNESSolve(user->snes_alg,NULL,X)); 224 } else if (event_list[i] == 1) { 225 /* Remove the fault */ 226 row_loc = 2*user->faultbus; col_loc = 2*user->faultbus+1; 227 val = -1/user->Rfault; 228 PetscCall(MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES)); 229 row_loc = 2*user->faultbus+1; col_loc = 2*user->faultbus; 230 val = -1/user->Rfault; 231 PetscCall(MatSetValues(user->Ybus,1,&row_loc,1,&col_loc,&val,ADD_VALUES)); 232 233 PetscCall(MatAssemblyBegin(user->Ybus,MAT_FINAL_ASSEMBLY)); 234 PetscCall(MatAssemblyEnd(user->Ybus,MAT_FINAL_ASSEMBLY)); 235 236 /* Solve the algebraic equations */ 237 PetscCall(SNESSolve(user->snes_alg,NULL,X)); 238 239 /* Check the VR derivatives and reset flags if needed */ 240 for (i=0; i < ngen; i++) { 241 Efd = xgen[idx+6]; 242 RF = xgen[idx+7]; 243 VR = xgen[idx+8]; 244 245 Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ 246 Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ 247 Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); 248 249 if (VRatmax[i]) { 250 fvalue = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; 251 if (fvalue < 0) { 252 VRatmax[i] = 0; 253 PetscCall(PetscPrintf(PETSC_COMM_SELF,"VR[%" PetscInt_FMT "]: dVR_dt went negative on fault clearing at time %g\n",i,(double)t)); 254 } 255 } 256 if (VRatmin[i]) { 257 fvalue = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; 258 259 if (fvalue > 0) { 260 VRatmin[i] = 0; 261 PetscCall(PetscPrintf(PETSC_COMM_SELF,"VR[%" PetscInt_FMT "]: dVR_dt went positive on fault clearing at time %g\n",i,(double)t)); 262 } 263 } 264 idx = idx+9; 265 } 266 } else { 267 idx = (event_list[i]-2)/2; 268 event_num = (event_list[i]-2)%2; 269 if (event_num == 0) { /* Max VR */ 270 if (!VRatmax[idx]) { 271 VRatmax[idx] = 1; 272 PetscCall(PetscPrintf(PETSC_COMM_SELF,"VR[%" PetscInt_FMT "]: hit upper limit at time %g\n",idx,(double)t)); 273 } 274 else { 275 VRatmax[idx] = 0; 276 PetscCall(PetscPrintf(PETSC_COMM_SELF,"VR[%" PetscInt_FMT "]: freeing variable as dVR_dt is negative at time %g\n",idx,(double)t)); 277 } 278 } else { 279 if (!VRatmin[idx]) { 280 VRatmin[idx] = 1; 281 PetscCall(PetscPrintf(PETSC_COMM_SELF,"VR[%" PetscInt_FMT "]: hit lower limit at time %g\n",idx,(double)t)); 282 } 283 else { 284 VRatmin[idx] = 0; 285 PetscCall(PetscPrintf(PETSC_COMM_SELF,"VR[%" PetscInt_FMT "]: freeing variable as dVR_dt is positive at time %g\n",idx,(double)t)); 286 } 287 } 288 } 289 } 290 291 PetscCall(VecRestoreArray(Xgen,&xgen)); 292 PetscCall(VecRestoreArray(Xnet,&xnet)); 293 294 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 295 296 PetscFunctionReturn(0); 297 } 298 299 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */ 300 PetscErrorCode dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar *Fr, PetscScalar *Fi) 301 { 302 PetscFunctionBegin; 303 *Fr = Fd*PetscSinScalar(delta) + Fq*PetscCosScalar(delta); 304 *Fi = -Fd*PetscCosScalar(delta) + Fq*PetscSinScalar(delta); 305 PetscFunctionReturn(0); 306 } 307 308 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */ 309 PetscErrorCode ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar *Fd, PetscScalar *Fq) 310 { 311 PetscFunctionBegin; 312 *Fd = Fr*PetscSinScalar(delta) - Fi*PetscCosScalar(delta); 313 *Fq = Fr*PetscCosScalar(delta) + Fi*PetscSinScalar(delta); 314 PetscFunctionReturn(0); 315 } 316 317 /* Saves the solution at each time to a matrix */ 318 PetscErrorCode SaveSolution(TS ts) 319 { 320 Userctx *user; 321 Vec X; 322 const PetscScalar *x; 323 PetscScalar *mat; 324 PetscInt idx; 325 PetscReal t; 326 327 PetscFunctionBegin; 328 PetscCall(TSGetApplicationContext(ts,&user)); 329 PetscCall(TSGetTime(ts,&t)); 330 PetscCall(TSGetSolution(ts,&X)); 331 idx = user->stepnum*(user->neqs_pgrid+1); 332 PetscCall(MatDenseGetArray(user->Sol,&mat)); 333 PetscCall(VecGetArrayRead(X,&x)); 334 mat[idx] = t; 335 PetscCall(PetscArraycpy(mat+idx+1,x,user->neqs_pgrid)); 336 PetscCall(MatDenseRestoreArray(user->Sol,&mat)); 337 PetscCall(VecRestoreArrayRead(X,&x)); 338 user->stepnum++; 339 PetscFunctionReturn(0); 340 } 341 342 PetscErrorCode SetInitialGuess(Vec X,Userctx *user) 343 { 344 Vec Xgen,Xnet; 345 PetscScalar *xgen; 346 const PetscScalar *xnet; 347 PetscInt i,idx=0; 348 PetscScalar Vr,Vi,IGr,IGi,Vm,Vm2; 349 PetscScalar Eqp,Edp,delta; 350 PetscScalar Efd,RF,VR; /* Exciter variables */ 351 PetscScalar Id,Iq; /* Generator dq axis currents */ 352 PetscScalar theta,Vd,Vq,SE; 353 354 PetscFunctionBegin; 355 M[0] = 2*H[0]/w_s; M[1] = 2*H[1]/w_s; M[2] = 2*H[2]/w_s; 356 D[0] = 0.1*M[0]; D[1] = 0.1*M[1]; D[2] = 0.1*M[2]; 357 358 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 359 360 /* Network subsystem initialization */ 361 PetscCall(VecCopy(user->V0,Xnet)); 362 363 /* Generator subsystem initialization */ 364 PetscCall(VecGetArrayWrite(Xgen,&xgen)); 365 PetscCall(VecGetArrayRead(Xnet,&xnet)); 366 367 for (i=0; i < ngen; i++) { 368 Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ 369 Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ 370 Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; 371 IGr = (Vr*PG[i] + Vi*QG[i])/Vm2; 372 IGi = (Vi*PG[i] - Vr*QG[i])/Vm2; 373 374 delta = PetscAtan2Real(Vi+Xq[i]*IGr,Vr-Xq[i]*IGi); /* Machine angle */ 375 376 theta = PETSC_PI/2.0 - delta; 377 378 Id = IGr*PetscCosScalar(theta) - IGi*PetscSinScalar(theta); /* d-axis stator current */ 379 Iq = IGr*PetscSinScalar(theta) + IGi*PetscCosScalar(theta); /* q-axis stator current */ 380 381 Vd = Vr*PetscCosScalar(theta) - Vi*PetscSinScalar(theta); 382 Vq = Vr*PetscSinScalar(theta) + Vi*PetscCosScalar(theta); 383 384 Edp = Vd + Rs[i]*Id - Xqp[i]*Iq; /* d-axis transient EMF */ 385 Eqp = Vq + Rs[i]*Iq + Xdp[i]*Id; /* q-axis transient EMF */ 386 387 TM[i] = PG[i]; 388 389 /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */ 390 xgen[idx] = Eqp; 391 xgen[idx+1] = Edp; 392 xgen[idx+2] = delta; 393 xgen[idx+3] = w_s; 394 395 idx = idx + 4; 396 397 xgen[idx] = Id; 398 xgen[idx+1] = Iq; 399 400 idx = idx + 2; 401 402 /* Exciter */ 403 Efd = Eqp + (Xd[i] - Xdp[i])*Id; 404 SE = k1[i]*PetscExpScalar(k2[i]*Efd); 405 VR = KE[i]*Efd + SE; 406 RF = KF[i]*Efd/TF[i]; 407 408 xgen[idx] = Efd; 409 xgen[idx+1] = RF; 410 xgen[idx+2] = VR; 411 412 Vref[i] = Vm + (VR/KA[i]); 413 414 VRatmin[i] = VRatmax[i] = 0; 415 416 idx = idx + 3; 417 } 418 PetscCall(VecRestoreArrayWrite(Xgen,&xgen)); 419 PetscCall(VecRestoreArrayRead(Xnet,&xnet)); 420 421 /* PetscCall(VecView(Xgen,0)); */ 422 PetscCall(DMCompositeGather(user->dmpgrid,INSERT_VALUES,X,Xgen,Xnet)); 423 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 424 PetscFunctionReturn(0); 425 } 426 427 /* Computes F = [f(x,y);g(x,y)] */ 428 PetscErrorCode ResidualFunction(Vec X, Vec F, Userctx *user) 429 { 430 Vec Xgen,Xnet,Fgen,Fnet; 431 const PetscScalar *xgen,*xnet; 432 PetscScalar *fgen,*fnet; 433 PetscInt i,idx=0; 434 PetscScalar Vr,Vi,Vm,Vm2; 435 PetscScalar Eqp,Edp,delta,w; /* Generator variables */ 436 PetscScalar Efd,RF,VR; /* Exciter variables */ 437 PetscScalar Id,Iq; /* Generator dq axis currents */ 438 PetscScalar Vd,Vq,SE; 439 PetscScalar IGr,IGi,IDr,IDi; 440 PetscScalar Zdq_inv[4],det; 441 PetscScalar PD,QD,Vm0,*v0; 442 PetscInt k; 443 444 PetscFunctionBegin; 445 PetscCall(VecZeroEntries(F)); 446 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 447 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Fgen,&Fnet)); 448 PetscCall(DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet)); 449 PetscCall(DMCompositeScatter(user->dmpgrid,F,Fgen,Fnet)); 450 451 /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here. 452 The generator current injection, IG, and load current injection, ID are added later 453 */ 454 /* Note that the values in Ybus are stored assuming the imaginary current balance 455 equation is ordered first followed by real current balance equation for each bus. 456 Thus imaginary current contribution goes in location 2*i, and 457 real current contribution in 2*i+1 458 */ 459 PetscCall(MatMult(user->Ybus,Xnet,Fnet)); 460 461 PetscCall(VecGetArrayRead(Xgen,&xgen)); 462 PetscCall(VecGetArrayRead(Xnet,&xnet)); 463 PetscCall(VecGetArrayWrite(Fgen,&fgen)); 464 PetscCall(VecGetArrayWrite(Fnet,&fnet)); 465 466 /* Generator subsystem */ 467 for (i=0; i < ngen; i++) { 468 Eqp = xgen[idx]; 469 Edp = xgen[idx+1]; 470 delta = xgen[idx+2]; 471 w = xgen[idx+3]; 472 Id = xgen[idx+4]; 473 Iq = xgen[idx+5]; 474 Efd = xgen[idx+6]; 475 RF = xgen[idx+7]; 476 VR = xgen[idx+8]; 477 478 /* Generator differential equations */ 479 fgen[idx] = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i]; 480 fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; 481 fgen[idx+2] = w - w_s; 482 fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i]; 483 484 Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ 485 Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ 486 487 PetscCall(ri2dq(Vr,Vi,delta,&Vd,&Vq)); 488 /* Algebraic equations for stator currents */ 489 det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i]; 490 491 Zdq_inv[0] = Rs[i]/det; 492 Zdq_inv[1] = Xqp[i]/det; 493 Zdq_inv[2] = -Xdp[i]/det; 494 Zdq_inv[3] = Rs[i]/det; 495 496 fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; 497 fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; 498 499 /* Add generator current injection to network */ 500 PetscCall(dq2ri(Id,Iq,delta,&IGr,&IGi)); 501 502 fnet[2*gbus[i]] -= IGi; 503 fnet[2*gbus[i]+1] -= IGr; 504 505 Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); 506 507 SE = k1[i]*PetscExpScalar(k2[i]*Efd); 508 509 /* Exciter differential equations */ 510 fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i]; 511 fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i]; 512 if (VRatmax[i]) fgen[idx+8] = VR - VRMAX[i]; 513 else if (VRatmin[i]) fgen[idx+8] = VRMIN[i] - VR; 514 else fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i]; 515 516 idx = idx + 9; 517 } 518 519 PetscCall(VecGetArray(user->V0,&v0)); 520 for (i=0; i < nload; i++) { 521 Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */ 522 Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */ 523 Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; 524 Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]); 525 PD = QD = 0.0; 526 for (k=0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]); 527 for (k=0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]); 528 529 /* Load currents */ 530 IDr = (PD*Vr + QD*Vi)/Vm2; 531 IDi = (-QD*Vr + PD*Vi)/Vm2; 532 533 fnet[2*lbus[i]] += IDi; 534 fnet[2*lbus[i]+1] += IDr; 535 } 536 PetscCall(VecRestoreArray(user->V0,&v0)); 537 538 PetscCall(VecRestoreArrayRead(Xgen,&xgen)); 539 PetscCall(VecRestoreArrayRead(Xnet,&xnet)); 540 PetscCall(VecRestoreArrayWrite(Fgen,&fgen)); 541 PetscCall(VecRestoreArrayWrite(Fnet,&fnet)); 542 543 PetscCall(DMCompositeGather(user->dmpgrid,INSERT_VALUES,F,Fgen,Fnet)); 544 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 545 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Fgen,&Fnet)); 546 PetscFunctionReturn(0); 547 } 548 549 /* f(x,y) 550 g(x,y) 551 */ 552 PetscErrorCode RHSFunction(TS ts,PetscReal t, Vec X, Vec F, void *ctx) 553 { 554 Userctx *user=(Userctx*)ctx; 555 556 PetscFunctionBegin; 557 user->t = t; 558 PetscCall(ResidualFunction(X,F,user)); 559 PetscFunctionReturn(0); 560 } 561 562 /* f(x,y) - \dot{x} 563 g(x,y) = 0 564 */ 565 PetscErrorCode IFunction(TS ts,PetscReal t, Vec X, Vec Xdot, Vec F, void *ctx) 566 { 567 PetscScalar *f; 568 const PetscScalar *xdot; 569 PetscInt i; 570 571 PetscFunctionBegin; 572 PetscCall(RHSFunction(ts,t,X,F,ctx)); 573 PetscCall(VecScale(F,-1.0)); 574 PetscCall(VecGetArray(F,&f)); 575 PetscCall(VecGetArrayRead(Xdot,&xdot)); 576 for (i=0;i < ngen;i++) { 577 f[9*i] += xdot[9*i]; 578 f[9*i+1] += xdot[9*i+1]; 579 f[9*i+2] += xdot[9*i+2]; 580 f[9*i+3] += xdot[9*i+3]; 581 f[9*i+6] += xdot[9*i+6]; 582 f[9*i+7] += xdot[9*i+7]; 583 f[9*i+8] += xdot[9*i+8]; 584 } 585 PetscCall(VecRestoreArray(F,&f)); 586 PetscCall(VecRestoreArrayRead(Xdot,&xdot)); 587 PetscFunctionReturn(0); 588 } 589 590 /* This function is used for solving the algebraic system only during fault on and 591 off times. It computes the entire F and then zeros out the part corresponding to 592 differential equations 593 F = [0;g(y)]; 594 */ 595 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx) 596 { 597 Userctx *user=(Userctx*)ctx; 598 PetscScalar *f; 599 PetscInt i; 600 601 PetscFunctionBegin; 602 PetscCall(ResidualFunction(X,F,user)); 603 PetscCall(VecGetArray(F,&f)); 604 for (i=0; i < ngen; i++) { 605 f[9*i] = 0; 606 f[9*i+1] = 0; 607 f[9*i+2] = 0; 608 f[9*i+3] = 0; 609 f[9*i+6] = 0; 610 f[9*i+7] = 0; 611 f[9*i+8] = 0; 612 } 613 PetscCall(VecRestoreArray(F,&f)); 614 PetscFunctionReturn(0); 615 } 616 617 PetscErrorCode PostStage(TS ts, PetscReal t, PetscInt i, Vec *X) 618 { 619 Userctx *user; 620 621 PetscFunctionBegin; 622 PetscCall(TSGetApplicationContext(ts,&user)); 623 PetscCall(SNESSolve(user->snes_alg,NULL,X[i])); 624 PetscFunctionReturn(0); 625 } 626 627 PetscErrorCode PostEvaluate(TS ts) 628 { 629 Userctx *user; 630 Vec X; 631 632 PetscFunctionBegin; 633 PetscCall(TSGetApplicationContext(ts,&user)); 634 PetscCall(TSGetSolution(ts,&X)); 635 PetscCall(SNESSolve(user->snes_alg,NULL,X)); 636 PetscFunctionReturn(0); 637 } 638 639 PetscErrorCode PreallocateJacobian(Mat J, Userctx *user) 640 { 641 PetscInt *d_nnz; 642 PetscInt i,idx=0,start=0; 643 PetscInt ncols; 644 645 PetscFunctionBegin; 646 PetscCall(PetscMalloc1(user->neqs_pgrid,&d_nnz)); 647 for (i=0; i<user->neqs_pgrid; i++) d_nnz[i] = 0; 648 /* Generator subsystem */ 649 for (i=0; i < ngen; i++) { 650 651 d_nnz[idx] += 3; 652 d_nnz[idx+1] += 2; 653 d_nnz[idx+2] += 2; 654 d_nnz[idx+3] += 5; 655 d_nnz[idx+4] += 6; 656 d_nnz[idx+5] += 6; 657 658 d_nnz[user->neqs_gen+2*gbus[i]] += 3; 659 d_nnz[user->neqs_gen+2*gbus[i]+1] += 3; 660 661 d_nnz[idx+6] += 2; 662 d_nnz[idx+7] += 2; 663 d_nnz[idx+8] += 5; 664 665 idx = idx + 9; 666 } 667 start = user->neqs_gen; 668 669 for (i=0; i < nbus; i++) { 670 PetscCall(MatGetRow(user->Ybus,2*i,&ncols,NULL,NULL)); 671 d_nnz[start+2*i] += ncols; 672 d_nnz[start+2*i+1] += ncols; 673 PetscCall(MatRestoreRow(user->Ybus,2*i,&ncols,NULL,NULL)); 674 } 675 PetscCall(MatSeqAIJSetPreallocation(J,0,d_nnz)); 676 PetscCall(PetscFree(d_nnz)); 677 PetscFunctionReturn(0); 678 } 679 680 /* 681 J = [df_dx, df_dy 682 dg_dx, dg_dy] 683 */ 684 PetscErrorCode ResidualJacobian(Vec X,Mat J,Mat B,void *ctx) 685 { 686 Userctx *user = (Userctx*)ctx; 687 Vec Xgen,Xnet; 688 const PetscScalar *xgen,*xnet; 689 PetscInt i,idx=0; 690 PetscScalar Vr,Vi,Vm,Vm2; 691 PetscScalar Eqp,Edp,delta; /* Generator variables */ 692 PetscScalar Efd; 693 PetscScalar Id,Iq; /* Generator dq axis currents */ 694 PetscScalar Vd,Vq; 695 PetscScalar val[10]; 696 PetscInt row[2],col[10]; 697 PetscInt net_start = user->neqs_gen; 698 PetscScalar Zdq_inv[4],det; 699 PetscScalar dVd_dVr,dVd_dVi,dVq_dVr,dVq_dVi,dVd_ddelta,dVq_ddelta; 700 PetscScalar dIGr_ddelta,dIGi_ddelta,dIGr_dId,dIGr_dIq,dIGi_dId,dIGi_dIq; 701 PetscScalar dSE_dEfd; 702 PetscScalar dVm_dVd,dVm_dVq,dVm_dVr,dVm_dVi; 703 PetscInt ncols; 704 const PetscInt *cols; 705 const PetscScalar *yvals; 706 PetscInt k; 707 PetscScalar PD,QD,Vm0,Vm4; 708 const PetscScalar *v0; 709 PetscScalar dPD_dVr,dPD_dVi,dQD_dVr,dQD_dVi; 710 PetscScalar dIDr_dVr,dIDr_dVi,dIDi_dVr,dIDi_dVi; 711 712 PetscFunctionBegin; 713 PetscCall(MatZeroEntries(B)); 714 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 715 PetscCall(DMCompositeScatter(user->dmpgrid,X,Xgen,Xnet)); 716 717 PetscCall(VecGetArrayRead(Xgen,&xgen)); 718 PetscCall(VecGetArrayRead(Xnet,&xnet)); 719 720 /* Generator subsystem */ 721 for (i=0; i < ngen; i++) { 722 Eqp = xgen[idx]; 723 Edp = xgen[idx+1]; 724 delta = xgen[idx+2]; 725 Id = xgen[idx+4]; 726 Iq = xgen[idx+5]; 727 Efd = xgen[idx+6]; 728 729 /* fgen[idx] = (-Eqp - (Xd[i] - Xdp[i])*Id + Efd)/Td0p[i]; */ 730 row[0] = idx; 731 col[0] = idx; col[1] = idx+4; col[2] = idx+6; 732 val[0] = -1/ Td0p[i]; val[1] = -(Xd[i] - Xdp[i])/ Td0p[i]; val[2] = 1/Td0p[i]; 733 734 PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES)); 735 736 /* fgen[idx+1] = (-Edp + (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */ 737 row[0] = idx + 1; 738 col[0] = idx + 1; col[1] = idx+5; 739 val[0] = -1/Tq0p[i]; val[1] = (Xq[i] - Xqp[i])/Tq0p[i]; 740 PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES)); 741 742 /* fgen[idx+2] = w - w_s; */ 743 row[0] = idx + 2; 744 col[0] = idx + 2; col[1] = idx + 3; 745 val[0] = 0; val[1] = 1; 746 PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES)); 747 748 /* fgen[idx+3] = (TM[i] - Edp*Id - Eqp*Iq - (Xqp[i] - Xdp[i])*Id*Iq - D[i]*(w - w_s))/M[i]; */ 749 row[0] = idx + 3; 750 col[0] = idx; col[1] = idx + 1; col[2] = idx + 3; col[3] = idx + 4; col[4] = idx + 5; 751 val[0] = -Iq/M[i]; val[1] = -Id/M[i]; val[2] = -D[i]/M[i]; val[3] = (-Edp - (Xqp[i]-Xdp[i])*Iq)/M[i]; val[4] = (-Eqp - (Xqp[i] - Xdp[i])*Id)/M[i]; 752 PetscCall(MatSetValues(J,1,row,5,col,val,INSERT_VALUES)); 753 754 Vr = xnet[2*gbus[i]]; /* Real part of generator terminal voltage */ 755 Vi = xnet[2*gbus[i]+1]; /* Imaginary part of the generator terminal voltage */ 756 PetscCall(ri2dq(Vr,Vi,delta,&Vd,&Vq)); 757 758 det = Rs[i]*Rs[i] + Xdp[i]*Xqp[i]; 759 760 Zdq_inv[0] = Rs[i]/det; 761 Zdq_inv[1] = Xqp[i]/det; 762 Zdq_inv[2] = -Xdp[i]/det; 763 Zdq_inv[3] = Rs[i]/det; 764 765 dVd_dVr = PetscSinScalar(delta); dVd_dVi = -PetscCosScalar(delta); 766 dVq_dVr = PetscCosScalar(delta); dVq_dVi = PetscSinScalar(delta); 767 dVd_ddelta = Vr*PetscCosScalar(delta) + Vi*PetscSinScalar(delta); 768 dVq_ddelta = -Vr*PetscSinScalar(delta) + Vi*PetscCosScalar(delta); 769 770 /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */ 771 row[0] = idx+4; 772 col[0] = idx; col[1] = idx+1; col[2] = idx + 2; 773 val[0] = -Zdq_inv[1]; val[1] = -Zdq_inv[0]; val[2] = Zdq_inv[0]*dVd_ddelta + Zdq_inv[1]*dVq_ddelta; 774 col[3] = idx + 4; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1; 775 val[3] = 1; val[4] = Zdq_inv[0]*dVd_dVr + Zdq_inv[1]*dVq_dVr; val[5] = Zdq_inv[0]*dVd_dVi + Zdq_inv[1]*dVq_dVi; 776 PetscCall(MatSetValues(J,1,row,6,col,val,INSERT_VALUES)); 777 778 /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */ 779 row[0] = idx+5; 780 col[0] = idx; col[1] = idx+1; col[2] = idx + 2; 781 val[0] = -Zdq_inv[3]; val[1] = -Zdq_inv[2]; val[2] = Zdq_inv[2]*dVd_ddelta + Zdq_inv[3]*dVq_ddelta; 782 col[3] = idx + 5; col[4] = net_start+2*gbus[i]; col[5] = net_start + 2*gbus[i]+1; 783 val[3] = 1; val[4] = Zdq_inv[2]*dVd_dVr + Zdq_inv[3]*dVq_dVr; val[5] = Zdq_inv[2]*dVd_dVi + Zdq_inv[3]*dVq_dVi; 784 PetscCall(MatSetValues(J,1,row,6,col,val,INSERT_VALUES)); 785 786 dIGr_ddelta = Id*PetscCosScalar(delta) - Iq*PetscSinScalar(delta); 787 dIGi_ddelta = Id*PetscSinScalar(delta) + Iq*PetscCosScalar(delta); 788 dIGr_dId = PetscSinScalar(delta); dIGr_dIq = PetscCosScalar(delta); 789 dIGi_dId = -PetscCosScalar(delta); dIGi_dIq = PetscSinScalar(delta); 790 791 /* fnet[2*gbus[i]] -= IGi; */ 792 row[0] = net_start + 2*gbus[i]; 793 col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5; 794 val[0] = -dIGi_ddelta; val[1] = -dIGi_dId; val[2] = -dIGi_dIq; 795 PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES)); 796 797 /* fnet[2*gbus[i]+1] -= IGr; */ 798 row[0] = net_start + 2*gbus[i]+1; 799 col[0] = idx+2; col[1] = idx + 4; col[2] = idx + 5; 800 val[0] = -dIGr_ddelta; val[1] = -dIGr_dId; val[2] = -dIGr_dIq; 801 PetscCall(MatSetValues(J,1,row,3,col,val,INSERT_VALUES)); 802 803 Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq); 804 805 /* fgen[idx+6] = (-KE[i]*Efd - SE + VR)/TE[i]; */ 806 /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */ 807 808 dSE_dEfd = k1[i]*k2[i]*PetscExpScalar(k2[i]*Efd); 809 810 row[0] = idx + 6; 811 col[0] = idx + 6; col[1] = idx + 8; 812 val[0] = (-KE[i] - dSE_dEfd)/TE[i]; val[1] = 1/TE[i]; 813 PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES)); 814 815 /* Exciter differential equations */ 816 817 /* fgen[idx+7] = (-RF + KF[i]*Efd/TF[i])/TF[i]; */ 818 row[0] = idx + 7; 819 col[0] = idx + 6; col[1] = idx + 7; 820 val[0] = (KF[i]/TF[i])/TF[i]; val[1] = -1/TF[i]; 821 PetscCall(MatSetValues(J,1,row,2,col,val,INSERT_VALUES)); 822 823 /* fgen[idx+8] = (-VR + KA[i]*RF - KA[i]*KF[i]*Efd/TF[i] + KA[i]*(Vref[i] - Vm))/TA[i]; */ 824 /* Vm = (Vd^2 + Vq^2)^0.5; */ 825 826 row[0] = idx + 8; 827 if (VRatmax[i]) { 828 col[0] = idx + 8; val[0] = 1.0; 829 PetscCall(MatSetValues(J,1,row,1,col,val,INSERT_VALUES)); 830 } else if (VRatmin[i]) { 831 col[0] = idx + 8; val[0] = -1.0; 832 PetscCall(MatSetValues(J,1,row,1,col,val,INSERT_VALUES)); 833 } else { 834 dVm_dVd = Vd/Vm; dVm_dVq = Vq/Vm; 835 dVm_dVr = dVm_dVd*dVd_dVr + dVm_dVq*dVq_dVr; 836 dVm_dVi = dVm_dVd*dVd_dVi + dVm_dVq*dVq_dVi; 837 row[0] = idx + 8; 838 col[0] = idx + 6; col[1] = idx + 7; col[2] = idx + 8; 839 val[0] = -(KA[i]*KF[i]/TF[i])/TA[i]; val[1] = KA[i]/TA[i]; val[2] = -1/TA[i]; 840 col[3] = net_start + 2*gbus[i]; col[4] = net_start + 2*gbus[i]+1; 841 val[3] = -KA[i]*dVm_dVr/TA[i]; val[4] = -KA[i]*dVm_dVi/TA[i]; 842 PetscCall(MatSetValues(J,1,row,5,col,val,INSERT_VALUES)); 843 } 844 idx = idx + 9; 845 } 846 847 for (i=0; i<nbus; i++) { 848 PetscCall(MatGetRow(user->Ybus,2*i,&ncols,&cols,&yvals)); 849 row[0] = net_start + 2*i; 850 for (k=0; k<ncols; k++) { 851 col[k] = net_start + cols[k]; 852 val[k] = yvals[k]; 853 } 854 PetscCall(MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES)); 855 PetscCall(MatRestoreRow(user->Ybus,2*i,&ncols,&cols,&yvals)); 856 857 PetscCall(MatGetRow(user->Ybus,2*i+1,&ncols,&cols,&yvals)); 858 row[0] = net_start + 2*i+1; 859 for (k=0; k<ncols; k++) { 860 col[k] = net_start + cols[k]; 861 val[k] = yvals[k]; 862 } 863 PetscCall(MatSetValues(J,1,row,ncols,col,val,INSERT_VALUES)); 864 PetscCall(MatRestoreRow(user->Ybus,2*i+1,&ncols,&cols,&yvals)); 865 } 866 867 PetscCall(MatAssemblyBegin(J,MAT_FLUSH_ASSEMBLY)); 868 PetscCall(MatAssemblyEnd(J,MAT_FLUSH_ASSEMBLY)); 869 870 PetscCall(VecGetArrayRead(user->V0,&v0)); 871 for (i=0; i < nload; i++) { 872 Vr = xnet[2*lbus[i]]; /* Real part of load bus voltage */ 873 Vi = xnet[2*lbus[i]+1]; /* Imaginary part of the load bus voltage */ 874 Vm = PetscSqrtScalar(Vr*Vr + Vi*Vi); Vm2 = Vm*Vm; Vm4 = Vm2*Vm2; 875 Vm0 = PetscSqrtScalar(v0[2*lbus[i]]*v0[2*lbus[i]] + v0[2*lbus[i]+1]*v0[2*lbus[i]+1]); 876 PD = QD = 0.0; 877 dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0; 878 for (k=0; k < ld_nsegsp[i]; k++) { 879 PD += ld_alphap[k]*PD0[i]*PetscPowScalar((Vm/Vm0),ld_betap[k]); 880 dPD_dVr += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vr*PetscPowScalar(Vm,(ld_betap[k]-2)); 881 dPD_dVi += ld_alphap[k]*ld_betap[k]*PD0[i]*PetscPowScalar((1/Vm0),ld_betap[k])*Vi*PetscPowScalar(Vm,(ld_betap[k]-2)); 882 } 883 for (k=0; k < ld_nsegsq[i]; k++) { 884 QD += ld_alphaq[k]*QD0[i]*PetscPowScalar((Vm/Vm0),ld_betaq[k]); 885 dQD_dVr += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vr*PetscPowScalar(Vm,(ld_betaq[k]-2)); 886 dQD_dVi += ld_alphaq[k]*ld_betaq[k]*QD0[i]*PetscPowScalar((1/Vm0),ld_betaq[k])*Vi*PetscPowScalar(Vm,(ld_betaq[k]-2)); 887 } 888 889 /* IDr = (PD*Vr + QD*Vi)/Vm2; */ 890 /* IDi = (-QD*Vr + PD*Vi)/Vm2; */ 891 892 dIDr_dVr = (dPD_dVr*Vr + dQD_dVr*Vi + PD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vr)/Vm4; 893 dIDr_dVi = (dPD_dVi*Vr + dQD_dVi*Vi + QD)/Vm2 - ((PD*Vr + QD*Vi)*2*Vi)/Vm4; 894 895 dIDi_dVr = (-dQD_dVr*Vr + dPD_dVr*Vi - QD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vr)/Vm4; 896 dIDi_dVi = (-dQD_dVi*Vr + dPD_dVi*Vi + PD)/Vm2 - ((-QD*Vr + PD*Vi)*2*Vi)/Vm4; 897 898 /* fnet[2*lbus[i]] += IDi; */ 899 row[0] = net_start + 2*lbus[i]; 900 col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1; 901 val[0] = dIDi_dVr; val[1] = dIDi_dVi; 902 PetscCall(MatSetValues(J,1,row,2,col,val,ADD_VALUES)); 903 /* fnet[2*lbus[i]+1] += IDr; */ 904 row[0] = net_start + 2*lbus[i]+1; 905 col[0] = net_start + 2*lbus[i]; col[1] = net_start + 2*lbus[i]+1; 906 val[0] = dIDr_dVr; val[1] = dIDr_dVi; 907 PetscCall(MatSetValues(J,1,row,2,col,val,ADD_VALUES)); 908 } 909 PetscCall(VecRestoreArrayRead(user->V0,&v0)); 910 911 PetscCall(VecRestoreArrayRead(Xgen,&xgen)); 912 PetscCall(VecRestoreArrayRead(Xnet,&xnet)); 913 914 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid,&Xgen,&Xnet)); 915 916 PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY)); 917 PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY)); 918 PetscFunctionReturn(0); 919 } 920 921 /* 922 J = [I, 0 923 dg_dx, dg_dy] 924 */ 925 PetscErrorCode AlgJacobian(SNES snes,Vec X,Mat A,Mat B,void *ctx) 926 { 927 Userctx *user=(Userctx*)ctx; 928 929 PetscFunctionBegin; 930 PetscCall(ResidualJacobian(X,A,B,ctx)); 931 PetscCall(MatSetOption(A,MAT_KEEP_NONZERO_PATTERN,PETSC_TRUE)); 932 PetscCall(MatZeroRowsIS(A,user->is_diff,1.0,NULL,NULL)); 933 PetscFunctionReturn(0); 934 } 935 936 /* 937 J = [-df_dx, -df_dy 938 dg_dx, dg_dy] 939 */ 940 941 PetscErrorCode RHSJacobian(TS ts,PetscReal t,Vec X,Mat A,Mat B,void *ctx) 942 { 943 Userctx *user=(Userctx*)ctx; 944 945 PetscFunctionBegin; 946 user->t = t; 947 948 PetscCall(ResidualJacobian(X,A,B,user)); 949 950 PetscFunctionReturn(0); 951 } 952 953 /* 954 J = [df_dx-aI, df_dy 955 dg_dx, dg_dy] 956 */ 957 958 PetscErrorCode IJacobian(TS ts,PetscReal t,Vec X,Vec Xdot,PetscReal a,Mat A,Mat B,Userctx *user) 959 { 960 PetscScalar atmp = (PetscScalar) a; 961 PetscInt i,row; 962 963 PetscFunctionBegin; 964 user->t = t; 965 966 PetscCall(RHSJacobian(ts,t,X,A,B,user)); 967 PetscCall(MatScale(B,-1.0)); 968 for (i=0;i < ngen;i++) { 969 row = 9*i; 970 PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES)); 971 row = 9*i+1; 972 PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES)); 973 row = 9*i+2; 974 PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES)); 975 row = 9*i+3; 976 PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES)); 977 row = 9*i+6; 978 PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES)); 979 row = 9*i+7; 980 PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES)); 981 row = 9*i+8; 982 PetscCall(MatSetValues(A,1,&row,1,&row,&atmp,ADD_VALUES)); 983 } 984 PetscCall(MatAssemblyBegin(A,MAT_FINAL_ASSEMBLY)); 985 PetscCall(MatAssemblyEnd(A,MAT_FINAL_ASSEMBLY)); 986 PetscFunctionReturn(0); 987 } 988 989 int main(int argc,char **argv) 990 { 991 TS ts; 992 SNES snes_alg; 993 PetscMPIInt size; 994 Userctx user; 995 PetscViewer Xview,Ybusview,viewer; 996 Vec X,F_alg; 997 Mat J,A; 998 PetscInt i,idx,*idx2; 999 Vec Xdot; 1000 PetscScalar *x,*mat,*amat; 1001 const PetscScalar *rmat; 1002 Vec vatol; 1003 PetscInt *direction; 1004 PetscBool *terminate; 1005 const PetscInt *idx3; 1006 PetscScalar *vatoli; 1007 PetscInt k; 1008 1009 PetscCall(PetscInitialize(&argc,&argv,"petscoptions",help)); 1010 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD,&size)); 1011 PetscCheck(size == 1,PETSC_COMM_WORLD,PETSC_ERR_WRONG_MPI_SIZE,"Only for sequential runs"); 1012 1013 user.neqs_gen = 9*ngen; /* # eqs. for generator subsystem */ 1014 user.neqs_net = 2*nbus; /* # eqs. for network subsystem */ 1015 user.neqs_pgrid = user.neqs_gen + user.neqs_net; 1016 1017 /* Create indices for differential and algebraic equations */ 1018 1019 PetscCall(PetscMalloc1(7*ngen,&idx2)); 1020 for (i=0; i<ngen; i++) { 1021 idx2[7*i] = 9*i; idx2[7*i+1] = 9*i+1; idx2[7*i+2] = 9*i+2; idx2[7*i+3] = 9*i+3; 1022 idx2[7*i+4] = 9*i+6; idx2[7*i+5] = 9*i+7; idx2[7*i+6] = 9*i+8; 1023 } 1024 PetscCall(ISCreateGeneral(PETSC_COMM_WORLD,7*ngen,idx2,PETSC_COPY_VALUES,&user.is_diff)); 1025 PetscCall(ISComplement(user.is_diff,0,user.neqs_pgrid,&user.is_alg)); 1026 PetscCall(PetscFree(idx2)); 1027 1028 /* Read initial voltage vector and Ybus */ 1029 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"X.bin",FILE_MODE_READ,&Xview)); 1030 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD,"Ybus.bin",FILE_MODE_READ,&Ybusview)); 1031 1032 PetscCall(VecCreate(PETSC_COMM_WORLD,&user.V0)); 1033 PetscCall(VecSetSizes(user.V0,PETSC_DECIDE,user.neqs_net)); 1034 PetscCall(VecLoad(user.V0,Xview)); 1035 1036 PetscCall(MatCreate(PETSC_COMM_WORLD,&user.Ybus)); 1037 PetscCall(MatSetSizes(user.Ybus,PETSC_DECIDE,PETSC_DECIDE,user.neqs_net,user.neqs_net)); 1038 PetscCall(MatSetType(user.Ybus,MATBAIJ)); 1039 /* PetscCall(MatSetBlockSize(user.Ybus,2)); */ 1040 PetscCall(MatLoad(user.Ybus,Ybusview)); 1041 1042 /* Set run time options */ 1043 PetscOptionsBegin(PETSC_COMM_WORLD,NULL,"Transient stability fault options",""); 1044 { 1045 user.tfaulton = 1.0; 1046 user.tfaultoff = 1.2; 1047 user.Rfault = 0.0001; 1048 user.setisdiff = PETSC_FALSE; 1049 user.semiexplicit = PETSC_FALSE; 1050 user.faultbus = 8; 1051 PetscCall(PetscOptionsReal("-tfaulton","","",user.tfaulton,&user.tfaulton,NULL)); 1052 PetscCall(PetscOptionsReal("-tfaultoff","","",user.tfaultoff,&user.tfaultoff,NULL)); 1053 PetscCall(PetscOptionsInt("-faultbus","","",user.faultbus,&user.faultbus,NULL)); 1054 user.t0 = 0.0; 1055 user.tmax = 5.0; 1056 PetscCall(PetscOptionsReal("-t0","","",user.t0,&user.t0,NULL)); 1057 PetscCall(PetscOptionsReal("-tmax","","",user.tmax,&user.tmax,NULL)); 1058 PetscCall(PetscOptionsBool("-setisdiff","","",user.setisdiff,&user.setisdiff,NULL)); 1059 PetscCall(PetscOptionsBool("-dae_semiexplicit","","",user.semiexplicit,&user.semiexplicit,NULL)); 1060 } 1061 PetscOptionsEnd(); 1062 1063 PetscCall(PetscViewerDestroy(&Xview)); 1064 PetscCall(PetscViewerDestroy(&Ybusview)); 1065 1066 /* Create DMs for generator and network subsystems */ 1067 PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_gen,1,1,NULL,&user.dmgen)); 1068 PetscCall(DMSetOptionsPrefix(user.dmgen,"dmgen_")); 1069 PetscCall(DMSetFromOptions(user.dmgen)); 1070 PetscCall(DMSetUp(user.dmgen)); 1071 PetscCall(DMDACreate1d(PETSC_COMM_WORLD,DM_BOUNDARY_NONE,user.neqs_net,1,1,NULL,&user.dmnet)); 1072 PetscCall(DMSetOptionsPrefix(user.dmnet,"dmnet_")); 1073 PetscCall(DMSetFromOptions(user.dmnet)); 1074 PetscCall(DMSetUp(user.dmnet)); 1075 /* Create a composite DM packer and add the two DMs */ 1076 PetscCall(DMCompositeCreate(PETSC_COMM_WORLD,&user.dmpgrid)); 1077 PetscCall(DMSetOptionsPrefix(user.dmpgrid,"pgrid_")); 1078 PetscCall(DMCompositeAddDM(user.dmpgrid,user.dmgen)); 1079 PetscCall(DMCompositeAddDM(user.dmpgrid,user.dmnet)); 1080 1081 PetscCall(DMCreateGlobalVector(user.dmpgrid,&X)); 1082 1083 PetscCall(MatCreate(PETSC_COMM_WORLD,&J)); 1084 PetscCall(MatSetSizes(J,PETSC_DECIDE,PETSC_DECIDE,user.neqs_pgrid,user.neqs_pgrid)); 1085 PetscCall(MatSetFromOptions(J)); 1086 PetscCall(PreallocateJacobian(J,&user)); 1087 1088 /* Create matrix to save solutions at each time step */ 1089 user.stepnum = 0; 1090 1091 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,1002,NULL,&user.Sol)); 1092 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1093 Create timestepping solver context 1094 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1095 PetscCall(TSCreate(PETSC_COMM_WORLD,&ts)); 1096 PetscCall(TSSetProblemType(ts,TS_NONLINEAR)); 1097 if (user.semiexplicit) { 1098 PetscCall(TSSetType(ts,TSRK)); 1099 PetscCall(TSSetRHSFunction(ts,NULL,RHSFunction,&user)); 1100 PetscCall(TSSetRHSJacobian(ts,J,J,RHSJacobian,&user)); 1101 } else { 1102 PetscCall(TSSetType(ts,TSCN)); 1103 PetscCall(TSSetEquationType(ts,TS_EQ_DAE_IMPLICIT_INDEX1)); 1104 PetscCall(TSSetIFunction(ts,NULL,(TSIFunction) IFunction,&user)); 1105 PetscCall(TSSetIJacobian(ts,J,J,(TSIJacobian)IJacobian,&user)); 1106 } 1107 PetscCall(TSSetApplicationContext(ts,&user)); 1108 1109 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 1110 Set initial conditions 1111 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 1112 PetscCall(SetInitialGuess(X,&user)); 1113 /* Just to set up the Jacobian structure */ 1114 1115 PetscCall(VecDuplicate(X,&Xdot)); 1116 PetscCall(IJacobian(ts,0.0,X,Xdot,0.0,J,J,&user)); 1117 PetscCall(VecDestroy(&Xdot)); 1118 1119 /* Save initial solution */ 1120 1121 idx=user.stepnum*(user.neqs_pgrid+1); 1122 PetscCall(MatDenseGetArray(user.Sol,&mat)); 1123 PetscCall(VecGetArray(X,&x)); 1124 1125 mat[idx] = 0.0; 1126 1127 PetscCall(PetscArraycpy(mat+idx+1,x,user.neqs_pgrid)); 1128 PetscCall(MatDenseRestoreArray(user.Sol,&mat)); 1129 PetscCall(VecRestoreArray(X,&x)); 1130 user.stepnum++; 1131 1132 PetscCall(TSSetMaxTime(ts,user.tmax)); 1133 PetscCall(TSSetExactFinalTime(ts,TS_EXACTFINALTIME_MATCHSTEP)); 1134 PetscCall(TSSetTimeStep(ts,0.01)); 1135 PetscCall(TSSetFromOptions(ts)); 1136 PetscCall(TSSetPostStep(ts,SaveSolution)); 1137 PetscCall(TSSetSolution(ts,X)); 1138 1139 PetscCall(PetscMalloc1((2*ngen+2),&direction)); 1140 PetscCall(PetscMalloc1((2*ngen+2),&terminate)); 1141 direction[0] = direction[1] = 1; 1142 terminate[0] = terminate[1] = PETSC_FALSE; 1143 for (i=0; i < ngen;i++) { 1144 direction[2+2*i] = -1; direction[2+2*i+1] = 1; 1145 terminate[2+2*i] = terminate[2+2*i+1] = PETSC_FALSE; 1146 } 1147 1148 PetscCall(TSSetEventHandler(ts,2*ngen+2,direction,terminate,EventFunction,PostEventFunction,(void*)&user)); 1149 1150 if (user.semiexplicit) { 1151 /* Use a semi-explicit approach with the time-stepping done by an explicit method and the 1152 algrebraic part solved via PostStage and PostEvaluate callbacks 1153 */ 1154 PetscCall(TSSetType(ts,TSRK)); 1155 PetscCall(TSSetPostStage(ts,PostStage)); 1156 PetscCall(TSSetPostEvaluate(ts,PostEvaluate)); 1157 } 1158 1159 if (user.setisdiff) { 1160 /* Create vector of absolute tolerances and set the algebraic part to infinity */ 1161 PetscCall(VecDuplicate(X,&vatol)); 1162 PetscCall(VecSet(vatol,100000.0)); 1163 PetscCall(VecGetArray(vatol,&vatoli)); 1164 PetscCall(ISGetIndices(user.is_diff,&idx3)); 1165 for (k=0; k < 7*ngen; k++) vatoli[idx3[k]] = 1e-2; 1166 PetscCall(VecRestoreArray(vatol,&vatoli)); 1167 } 1168 1169 /* Create the nonlinear solver for solving the algebraic system */ 1170 /* Note that although the algebraic system needs to be solved only for 1171 Idq and V, we reuse the entire system including xgen. The xgen 1172 variables are held constant by setting their residuals to 0 and 1173 putting a 1 on the Jacobian diagonal for xgen rows 1174 */ 1175 1176 PetscCall(VecDuplicate(X,&F_alg)); 1177 PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes_alg)); 1178 PetscCall(SNESSetFunction(snes_alg,F_alg,AlgFunction,&user)); 1179 PetscCall(SNESSetJacobian(snes_alg,J,J,AlgJacobian,&user)); 1180 PetscCall(SNESSetFromOptions(snes_alg)); 1181 1182 user.snes_alg=snes_alg; 1183 1184 /* Solve */ 1185 PetscCall(TSSolve(ts,X)); 1186 1187 PetscCall(MatAssemblyBegin(user.Sol,MAT_FINAL_ASSEMBLY)); 1188 PetscCall(MatAssemblyEnd(user.Sol,MAT_FINAL_ASSEMBLY)); 1189 1190 PetscCall(MatCreateSeqDense(PETSC_COMM_SELF,user.neqs_pgrid+1,user.stepnum,NULL,&A)); 1191 PetscCall(MatDenseGetArrayRead(user.Sol,&rmat)); 1192 PetscCall(MatDenseGetArray(A,&amat)); 1193 PetscCall(PetscArraycpy(amat,rmat,user.stepnum*(user.neqs_pgrid+1))); 1194 PetscCall(MatDenseRestoreArray(A,&amat)); 1195 PetscCall(MatDenseRestoreArrayRead(user.Sol,&rmat)); 1196 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF,"out.bin",FILE_MODE_WRITE,&viewer)); 1197 PetscCall(MatView(A,viewer)); 1198 PetscCall(PetscViewerDestroy(&viewer)); 1199 PetscCall(MatDestroy(&A)); 1200 1201 PetscCall(PetscFree(direction)); 1202 PetscCall(PetscFree(terminate)); 1203 PetscCall(SNESDestroy(&snes_alg)); 1204 PetscCall(VecDestroy(&F_alg)); 1205 PetscCall(MatDestroy(&J)); 1206 PetscCall(MatDestroy(&user.Ybus)); 1207 PetscCall(MatDestroy(&user.Sol)); 1208 PetscCall(VecDestroy(&X)); 1209 PetscCall(VecDestroy(&user.V0)); 1210 PetscCall(DMDestroy(&user.dmgen)); 1211 PetscCall(DMDestroy(&user.dmnet)); 1212 PetscCall(DMDestroy(&user.dmpgrid)); 1213 PetscCall(ISDestroy(&user.is_diff)); 1214 PetscCall(ISDestroy(&user.is_alg)); 1215 PetscCall(TSDestroy(&ts)); 1216 if (user.setisdiff) { 1217 PetscCall(VecDestroy(&vatol)); 1218 } 1219 PetscCall(PetscFinalize()); 1220 return 0; 1221 } 1222 1223 /*TEST 1224 1225 build: 1226 requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 1227 1228 test: 1229 suffix: implicit 1230 args: -ts_monitor -snes_monitor_short 1231 localrunfiles: petscoptions X.bin Ybus.bin 1232 1233 test: 1234 suffix: semiexplicit 1235 args: -ts_monitor -snes_monitor_short -dae_semiexplicit -ts_rk_type 2a 1236 localrunfiles: petscoptions X.bin Ybus.bin 1237 1238 test: 1239 suffix: steprestart 1240 # needs ARKIMEX methods with all implicit stages since the mass matrix is not the identity 1241 args: -ts_monitor -snes_monitor_short -ts_type arkimex -ts_arkimex_type prssp2 1242 localrunfiles: petscoptions X.bin Ybus.bin 1243 1244 TEST*/ 1245