1 2 static char help[] = "Sensitivity analysis applied in 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. See ex9bus.c for details. 11 The system has 'jumps' due to faults, thus the time interval is split into multiple sections, and TSSolve is called for each of them. But TSAdjointSolve only needs to be called once since the whole trajectory has been saved in the forward run. 12 The code computes the sensitivity of a final state w.r.t. initial conditions. 13 */ 14 15 #include <petscts.h> 16 #include <petscdm.h> 17 #include <petscdmda.h> 18 #include <petscdmcomposite.h> 19 20 #define freq 60 21 #define w_s (2 * PETSC_PI * freq) 22 23 /* Sizes and indices */ 24 const PetscInt nbus = 9; /* Number of network buses */ 25 const PetscInt ngen = 3; /* Number of generators */ 26 const PetscInt nload = 3; /* Number of loads */ 27 const PetscInt gbus[3] = {0, 1, 2}; /* Buses at which generators are incident */ 28 const PetscInt lbus[3] = {4, 5, 7}; /* Buses at which loads are incident */ 29 30 /* Generator real and reactive powers (found via loadflow) */ 31 const PetscScalar PG[3] = {0.716786142395021, 1.630000000000000, 0.850000000000000}; 32 const PetscScalar QG[3] = {0.270702180178785, 0.066120127797275, -0.108402221791588}; 33 /* Generator constants */ 34 const PetscScalar H[3] = {23.64, 6.4, 3.01}; /* Inertia constant */ 35 const PetscScalar Rs[3] = {0.0, 0.0, 0.0}; /* Stator Resistance */ 36 const PetscScalar Xd[3] = {0.146, 0.8958, 1.3125}; /* d-axis reactance */ 37 const PetscScalar Xdp[3] = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */ 38 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 */ 39 const PetscScalar Xqp[3] = {0.0969, 0.1969, 0.25}; /* q-axis transient reactance */ 40 const PetscScalar Td0p[3] = {8.96, 6.0, 5.89}; /* d-axis open circuit time constant */ 41 const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6}; /* q-axis open circuit time constant */ 42 PetscScalar M[3]; /* M = 2*H/w_s */ 43 PetscScalar D[3]; /* D = 0.1*M */ 44 45 PetscScalar TM[3]; /* Mechanical Torque */ 46 /* Exciter system constants */ 47 const PetscScalar KA[3] = {20.0, 20.0, 20.0}; /* Voltage regulartor gain constant */ 48 const PetscScalar TA[3] = {0.2, 0.2, 0.2}; /* Voltage regulator time constant */ 49 const PetscScalar KE[3] = {1.0, 1.0, 1.0}; /* Exciter gain constant */ 50 const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */ 51 const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */ 52 const PetscScalar TF[3] = {0.35, 0.35, 0.35}; /* Feedback stabilizer time constant */ 53 const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039}; 54 const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */ 55 56 PetscScalar Vref[3]; 57 /* Load constants 58 We use a composite load model that describes the load and reactive powers at each time instant as follows 59 P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i 60 Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i 61 where 62 ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads 63 ld_alphap,ld_alphap - Percentage contribution (weights) or loads 64 P_D0 - Real power load 65 Q_D0 - Reactive power load 66 V_m(t) - Voltage magnitude at time t 67 V_m0 - Voltage magnitude at t = 0 68 ld_betap, ld_betaq - exponents describing the load model for real and reactive part 69 70 Note: All loads have the same characteristic currently. 71 */ 72 const PetscScalar PD0[3] = {1.25, 0.9, 1.0}; 73 const PetscScalar QD0[3] = {0.5, 0.3, 0.35}; 74 const PetscInt ld_nsegsp[3] = {3, 3, 3}; 75 const PetscScalar ld_alphap[3] = {1.0, 0.0, 0.0}; 76 const PetscScalar ld_betap[3] = {2.0, 1.0, 0.0}; 77 const PetscInt ld_nsegsq[3] = {3, 3, 3}; 78 const PetscScalar ld_alphaq[3] = {1.0, 0.0, 0.0}; 79 const PetscScalar ld_betaq[3] = {2.0, 1.0, 0.0}; 80 81 typedef struct { 82 DM dmgen, dmnet; /* DMs to manage generator and network subsystem */ 83 DM dmpgrid; /* Composite DM to manage the entire power grid */ 84 Mat Ybus; /* Network admittance matrix */ 85 Vec V0; /* Initial voltage vector (Power flow solution) */ 86 PetscReal tfaulton, tfaultoff; /* Fault on and off times */ 87 PetscInt faultbus; /* Fault bus */ 88 PetscScalar Rfault; 89 PetscReal t0, tmax; 90 PetscInt neqs_gen, neqs_net, neqs_pgrid; 91 PetscBool alg_flg; 92 PetscReal t; 93 IS is_diff; /* indices for differential equations */ 94 IS is_alg; /* indices for algebraic equations */ 95 } Userctx; 96 97 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */ 98 PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi) 99 { 100 PetscFunctionBegin; 101 *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta); 102 *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta); 103 PetscFunctionReturn(PETSC_SUCCESS); 104 } 105 106 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */ 107 PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq) 108 { 109 PetscFunctionBegin; 110 *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta); 111 *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta); 112 PetscFunctionReturn(PETSC_SUCCESS); 113 } 114 115 PetscErrorCode SetInitialGuess(Vec X, Userctx *user) 116 { 117 Vec Xgen, Xnet; 118 PetscScalar *xgen, *xnet; 119 PetscInt i, idx = 0; 120 PetscScalar Vr, Vi, IGr, IGi, Vm, Vm2; 121 PetscScalar Eqp, Edp, delta; 122 PetscScalar Efd, RF, VR; /* Exciter variables */ 123 PetscScalar Id, Iq; /* Generator dq axis currents */ 124 PetscScalar theta, Vd, Vq, SE; 125 126 PetscFunctionBegin; 127 M[0] = 2 * H[0] / w_s; 128 M[1] = 2 * H[1] / w_s; 129 M[2] = 2 * H[2] / w_s; 130 D[0] = 0.1 * M[0]; 131 D[1] = 0.1 * M[1]; 132 D[2] = 0.1 * M[2]; 133 134 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet)); 135 136 /* Network subsystem initialization */ 137 PetscCall(VecCopy(user->V0, Xnet)); 138 139 /* Generator subsystem initialization */ 140 PetscCall(VecGetArray(Xgen, &xgen)); 141 PetscCall(VecGetArray(Xnet, &xnet)); 142 143 for (i = 0; i < ngen; i++) { 144 Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */ 145 Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */ 146 Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi); 147 Vm2 = Vm * Vm; 148 IGr = (Vr * PG[i] + Vi * QG[i]) / Vm2; 149 IGi = (Vi * PG[i] - Vr * QG[i]) / Vm2; 150 151 delta = PetscAtan2Real(Vi + Xq[i] * IGr, Vr - Xq[i] * IGi); /* Machine angle */ 152 153 theta = PETSC_PI / 2.0 - delta; 154 155 Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta); /* d-axis stator current */ 156 Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta); /* q-axis stator current */ 157 158 Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta); 159 Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta); 160 161 Edp = Vd + Rs[i] * Id - Xqp[i] * Iq; /* d-axis transient EMF */ 162 Eqp = Vq + Rs[i] * Iq + Xdp[i] * Id; /* q-axis transient EMF */ 163 164 TM[i] = PG[i]; 165 166 /* The generator variables are ordered as [Eqp,Edp,delta,w,Id,Iq] */ 167 xgen[idx] = Eqp; 168 xgen[idx + 1] = Edp; 169 xgen[idx + 2] = delta; 170 xgen[idx + 3] = w_s; 171 172 idx = idx + 4; 173 174 xgen[idx] = Id; 175 xgen[idx + 1] = Iq; 176 177 idx = idx + 2; 178 179 /* Exciter */ 180 Efd = Eqp + (Xd[i] - Xdp[i]) * Id; 181 SE = k1[i] * PetscExpScalar(k2[i] * Efd); 182 VR = KE[i] * Efd + SE; 183 RF = KF[i] * Efd / TF[i]; 184 185 xgen[idx] = Efd; 186 xgen[idx + 1] = RF; 187 xgen[idx + 2] = VR; 188 189 Vref[i] = Vm + (VR / KA[i]); 190 191 idx = idx + 3; 192 } 193 194 PetscCall(VecRestoreArray(Xgen, &xgen)); 195 PetscCall(VecRestoreArray(Xnet, &xnet)); 196 197 /* PetscCall(VecView(Xgen,0)); */ 198 PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, X, Xgen, Xnet)); 199 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet)); 200 PetscFunctionReturn(PETSC_SUCCESS); 201 } 202 203 /* Computes F = [f(x,y);g(x,y)] */ 204 PetscErrorCode ResidualFunction(SNES snes, Vec X, Vec F, Userctx *user) 205 { 206 Vec Xgen, Xnet, Fgen, Fnet; 207 PetscScalar *xgen, *xnet, *fgen, *fnet; 208 PetscInt i, idx = 0; 209 PetscScalar Vr, Vi, Vm, Vm2; 210 PetscScalar Eqp, Edp, delta, w; /* Generator variables */ 211 PetscScalar Efd, RF, VR; /* Exciter variables */ 212 PetscScalar Id, Iq; /* Generator dq axis currents */ 213 PetscScalar Vd, Vq, SE; 214 PetscScalar IGr, IGi, IDr, IDi; 215 PetscScalar Zdq_inv[4], det; 216 PetscScalar PD, QD, Vm0, *v0; 217 PetscInt k; 218 219 PetscFunctionBegin; 220 PetscCall(VecZeroEntries(F)); 221 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet)); 222 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Fgen, &Fnet)); 223 PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet)); 224 PetscCall(DMCompositeScatter(user->dmpgrid, F, Fgen, Fnet)); 225 226 /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here. 227 The generator current injection, IG, and load current injection, ID are added later 228 */ 229 /* Note that the values in Ybus are stored assuming the imaginary current balance 230 equation is ordered first followed by real current balance equation for each bus. 231 Thus imaginary current contribution goes in location 2*i, and 232 real current contribution in 2*i+1 233 */ 234 PetscCall(MatMult(user->Ybus, Xnet, Fnet)); 235 236 PetscCall(VecGetArray(Xgen, &xgen)); 237 PetscCall(VecGetArray(Xnet, &xnet)); 238 PetscCall(VecGetArray(Fgen, &fgen)); 239 PetscCall(VecGetArray(Fnet, &fnet)); 240 241 /* Generator subsystem */ 242 for (i = 0; i < ngen; i++) { 243 Eqp = xgen[idx]; 244 Edp = xgen[idx + 1]; 245 delta = xgen[idx + 2]; 246 w = xgen[idx + 3]; 247 Id = xgen[idx + 4]; 248 Iq = xgen[idx + 5]; 249 Efd = xgen[idx + 6]; 250 RF = xgen[idx + 7]; 251 VR = xgen[idx + 8]; 252 253 /* Generator differential equations */ 254 fgen[idx] = (Eqp + (Xd[i] - Xdp[i]) * Id - Efd) / Td0p[i]; 255 fgen[idx + 1] = (Edp - (Xq[i] - Xqp[i]) * Iq) / Tq0p[i]; 256 fgen[idx + 2] = -w + w_s; 257 fgen[idx + 3] = (-TM[i] + Edp * Id + Eqp * Iq + (Xqp[i] - Xdp[i]) * Id * Iq + D[i] * (w - w_s)) / M[i]; 258 259 Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */ 260 Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */ 261 262 PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq)); 263 /* Algebraic equations for stator currents */ 264 265 det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i]; 266 267 Zdq_inv[0] = Rs[i] / det; 268 Zdq_inv[1] = Xqp[i] / det; 269 Zdq_inv[2] = -Xdp[i] / det; 270 Zdq_inv[3] = Rs[i] / det; 271 272 fgen[idx + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id; 273 fgen[idx + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq; 274 275 /* Add generator current injection to network */ 276 PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi)); 277 278 fnet[2 * gbus[i]] -= IGi; 279 fnet[2 * gbus[i] + 1] -= IGr; 280 281 Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq); 282 283 SE = k1[i] * PetscExpScalar(k2[i] * Efd); 284 285 /* Exciter differential equations */ 286 fgen[idx + 6] = (KE[i] * Efd + SE - VR) / TE[i]; 287 fgen[idx + 7] = (RF - KF[i] * Efd / TF[i]) / TF[i]; 288 fgen[idx + 8] = (VR - KA[i] * RF + KA[i] * KF[i] * Efd / TF[i] - KA[i] * (Vref[i] - Vm)) / TA[i]; 289 290 idx = idx + 9; 291 } 292 293 PetscCall(VecGetArray(user->V0, &v0)); 294 for (i = 0; i < nload; i++) { 295 Vr = xnet[2 * lbus[i]]; /* Real part of load bus voltage */ 296 Vi = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */ 297 Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi); 298 Vm2 = Vm * Vm; 299 Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]); 300 PD = QD = 0.0; 301 for (k = 0; k < ld_nsegsp[i]; k++) PD += ld_alphap[k] * PD0[i] * PetscPowScalar((Vm / Vm0), ld_betap[k]); 302 for (k = 0; k < ld_nsegsq[i]; k++) QD += ld_alphaq[k] * QD0[i] * PetscPowScalar((Vm / Vm0), ld_betaq[k]); 303 304 /* Load currents */ 305 IDr = (PD * Vr + QD * Vi) / Vm2; 306 IDi = (-QD * Vr + PD * Vi) / Vm2; 307 308 fnet[2 * lbus[i]] += IDi; 309 fnet[2 * lbus[i] + 1] += IDr; 310 } 311 PetscCall(VecRestoreArray(user->V0, &v0)); 312 313 PetscCall(VecRestoreArray(Xgen, &xgen)); 314 PetscCall(VecRestoreArray(Xnet, &xnet)); 315 PetscCall(VecRestoreArray(Fgen, &fgen)); 316 PetscCall(VecRestoreArray(Fnet, &fnet)); 317 318 PetscCall(DMCompositeGather(user->dmpgrid, INSERT_VALUES, F, Fgen, Fnet)); 319 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet)); 320 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Fgen, &Fnet)); 321 PetscFunctionReturn(PETSC_SUCCESS); 322 } 323 324 /* \dot{x} - f(x,y) 325 g(x,y) = 0 326 */ 327 PetscErrorCode IFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user) 328 { 329 SNES snes; 330 PetscScalar *f; 331 const PetscScalar *xdot; 332 PetscInt i; 333 334 PetscFunctionBegin; 335 user->t = t; 336 337 PetscCall(TSGetSNES(ts, &snes)); 338 PetscCall(ResidualFunction(snes, X, F, user)); 339 PetscCall(VecGetArray(F, &f)); 340 PetscCall(VecGetArrayRead(Xdot, &xdot)); 341 for (i = 0; i < ngen; i++) { 342 f[9 * i] += xdot[9 * i]; 343 f[9 * i + 1] += xdot[9 * i + 1]; 344 f[9 * i + 2] += xdot[9 * i + 2]; 345 f[9 * i + 3] += xdot[9 * i + 3]; 346 f[9 * i + 6] += xdot[9 * i + 6]; 347 f[9 * i + 7] += xdot[9 * i + 7]; 348 f[9 * i + 8] += xdot[9 * i + 8]; 349 } 350 PetscCall(VecRestoreArray(F, &f)); 351 PetscCall(VecRestoreArrayRead(Xdot, &xdot)); 352 PetscFunctionReturn(PETSC_SUCCESS); 353 } 354 355 /* This function is used for solving the algebraic system only during fault on and 356 off times. It computes the entire F and then zeros out the part corresponding to 357 differential equations 358 F = [0;g(y)]; 359 */ 360 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx) 361 { 362 Userctx *user = (Userctx *)ctx; 363 PetscScalar *f; 364 PetscInt i; 365 366 PetscFunctionBegin; 367 PetscCall(ResidualFunction(snes, X, F, user)); 368 PetscCall(VecGetArray(F, &f)); 369 for (i = 0; i < ngen; i++) { 370 f[9 * i] = 0; 371 f[9 * i + 1] = 0; 372 f[9 * i + 2] = 0; 373 f[9 * i + 3] = 0; 374 f[9 * i + 6] = 0; 375 f[9 * i + 7] = 0; 376 f[9 * i + 8] = 0; 377 } 378 PetscCall(VecRestoreArray(F, &f)); 379 PetscFunctionReturn(PETSC_SUCCESS); 380 } 381 382 PetscErrorCode PreallocateJacobian(Mat J, Userctx *user) 383 { 384 PetscInt *d_nnz; 385 PetscInt i, idx = 0, start = 0; 386 PetscInt ncols; 387 388 PetscFunctionBegin; 389 PetscCall(PetscMalloc1(user->neqs_pgrid, &d_nnz)); 390 for (i = 0; i < user->neqs_pgrid; i++) d_nnz[i] = 0; 391 /* Generator subsystem */ 392 for (i = 0; i < ngen; i++) { 393 d_nnz[idx] += 3; 394 d_nnz[idx + 1] += 2; 395 d_nnz[idx + 2] += 2; 396 d_nnz[idx + 3] += 5; 397 d_nnz[idx + 4] += 6; 398 d_nnz[idx + 5] += 6; 399 400 d_nnz[user->neqs_gen + 2 * gbus[i]] += 3; 401 d_nnz[user->neqs_gen + 2 * gbus[i] + 1] += 3; 402 403 d_nnz[idx + 6] += 2; 404 d_nnz[idx + 7] += 2; 405 d_nnz[idx + 8] += 5; 406 407 idx = idx + 9; 408 } 409 410 start = user->neqs_gen; 411 412 for (i = 0; i < nbus; i++) { 413 PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, NULL, NULL)); 414 d_nnz[start + 2 * i] += ncols; 415 d_nnz[start + 2 * i + 1] += ncols; 416 PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, NULL, NULL)); 417 } 418 419 PetscCall(MatSeqAIJSetPreallocation(J, 0, d_nnz)); 420 421 PetscCall(PetscFree(d_nnz)); 422 PetscFunctionReturn(PETSC_SUCCESS); 423 } 424 425 /* 426 J = [-df_dx, -df_dy 427 dg_dx, dg_dy] 428 */ 429 PetscErrorCode ResidualJacobian(SNES snes, Vec X, Mat J, Mat B, void *ctx) 430 { 431 Userctx *user = (Userctx *)ctx; 432 Vec Xgen, Xnet; 433 PetscScalar *xgen, *xnet; 434 PetscInt i, idx = 0; 435 PetscScalar Vr, Vi, Vm, Vm2; 436 PetscScalar Eqp, Edp, delta; /* Generator variables */ 437 PetscScalar Efd; /* Exciter variables */ 438 PetscScalar Id, Iq; /* Generator dq axis currents */ 439 PetscScalar Vd, Vq; 440 PetscScalar val[10]; 441 PetscInt row[2], col[10]; 442 PetscInt net_start = user->neqs_gen; 443 PetscScalar Zdq_inv[4], det; 444 PetscScalar dVd_dVr, dVd_dVi, dVq_dVr, dVq_dVi, dVd_ddelta, dVq_ddelta; 445 PetscScalar dIGr_ddelta, dIGi_ddelta, dIGr_dId, dIGr_dIq, dIGi_dId, dIGi_dIq; 446 PetscScalar dSE_dEfd; 447 PetscScalar dVm_dVd, dVm_dVq, dVm_dVr, dVm_dVi; 448 PetscInt ncols; 449 const PetscInt *cols; 450 const PetscScalar *yvals; 451 PetscInt k; 452 PetscScalar PD, QD, Vm0, *v0, Vm4; 453 PetscScalar dPD_dVr, dPD_dVi, dQD_dVr, dQD_dVi; 454 PetscScalar dIDr_dVr, dIDr_dVi, dIDi_dVr, dIDi_dVi; 455 456 PetscFunctionBegin; 457 PetscCall(MatZeroEntries(B)); 458 PetscCall(DMCompositeGetLocalVectors(user->dmpgrid, &Xgen, &Xnet)); 459 PetscCall(DMCompositeScatter(user->dmpgrid, X, Xgen, Xnet)); 460 461 PetscCall(VecGetArray(Xgen, &xgen)); 462 PetscCall(VecGetArray(Xnet, &xnet)); 463 464 /* Generator subsystem */ 465 for (i = 0; i < ngen; i++) { 466 Eqp = xgen[idx]; 467 Edp = xgen[idx + 1]; 468 delta = xgen[idx + 2]; 469 Id = xgen[idx + 4]; 470 Iq = xgen[idx + 5]; 471 Efd = xgen[idx + 6]; 472 473 /* fgen[idx] = (Eqp + (Xd[i] - Xdp[i])*Id - Efd)/Td0p[i]; */ 474 row[0] = idx; 475 col[0] = idx; 476 col[1] = idx + 4; 477 col[2] = idx + 6; 478 val[0] = 1 / Td0p[i]; 479 val[1] = (Xd[i] - Xdp[i]) / Td0p[i]; 480 val[2] = -1 / Td0p[i]; 481 482 PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES)); 483 484 /* fgen[idx+1] = (Edp - (Xq[i] - Xqp[i])*Iq)/Tq0p[i]; */ 485 row[0] = idx + 1; 486 col[0] = idx + 1; 487 col[1] = idx + 5; 488 val[0] = 1 / Tq0p[i]; 489 val[1] = -(Xq[i] - Xqp[i]) / Tq0p[i]; 490 PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES)); 491 492 /* fgen[idx+2] = - w + w_s; */ 493 row[0] = idx + 2; 494 col[0] = idx + 2; 495 col[1] = idx + 3; 496 val[0] = 0; 497 val[1] = -1; 498 PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES)); 499 500 /* fgen[idx+3] = (-TM[i] + Edp*Id + Eqp*Iq + (Xqp[i] - Xdp[i])*Id*Iq + D[i]*(w - w_s))/M[i]; */ 501 row[0] = idx + 3; 502 col[0] = idx; 503 col[1] = idx + 1; 504 col[2] = idx + 3; 505 col[3] = idx + 4; 506 col[4] = idx + 5; 507 val[0] = Iq / M[i]; 508 val[1] = Id / M[i]; 509 val[2] = D[i] / M[i]; 510 val[3] = (Edp + (Xqp[i] - Xdp[i]) * Iq) / M[i]; 511 val[4] = (Eqp + (Xqp[i] - Xdp[i]) * Id) / M[i]; 512 PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES)); 513 514 Vr = xnet[2 * gbus[i]]; /* Real part of generator terminal voltage */ 515 Vi = xnet[2 * gbus[i] + 1]; /* Imaginary part of the generator terminal voltage */ 516 PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq)); 517 518 det = Rs[i] * Rs[i] + Xdp[i] * Xqp[i]; 519 520 Zdq_inv[0] = Rs[i] / det; 521 Zdq_inv[1] = Xqp[i] / det; 522 Zdq_inv[2] = -Xdp[i] / det; 523 Zdq_inv[3] = Rs[i] / det; 524 525 dVd_dVr = PetscSinScalar(delta); 526 dVd_dVi = -PetscCosScalar(delta); 527 dVq_dVr = PetscCosScalar(delta); 528 dVq_dVi = PetscSinScalar(delta); 529 dVd_ddelta = Vr * PetscCosScalar(delta) + Vi * PetscSinScalar(delta); 530 dVq_ddelta = -Vr * PetscSinScalar(delta) + Vi * PetscCosScalar(delta); 531 532 /* fgen[idx+4] = Zdq_inv[0]*(-Edp + Vd) + Zdq_inv[1]*(-Eqp + Vq) + Id; */ 533 row[0] = idx + 4; 534 col[0] = idx; 535 col[1] = idx + 1; 536 col[2] = idx + 2; 537 val[0] = -Zdq_inv[1]; 538 val[1] = -Zdq_inv[0]; 539 val[2] = Zdq_inv[0] * dVd_ddelta + Zdq_inv[1] * dVq_ddelta; 540 col[3] = idx + 4; 541 col[4] = net_start + 2 * gbus[i]; 542 col[5] = net_start + 2 * gbus[i] + 1; 543 val[3] = 1; 544 val[4] = Zdq_inv[0] * dVd_dVr + Zdq_inv[1] * dVq_dVr; 545 val[5] = Zdq_inv[0] * dVd_dVi + Zdq_inv[1] * dVq_dVi; 546 PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES)); 547 548 /* fgen[idx+5] = Zdq_inv[2]*(-Edp + Vd) + Zdq_inv[3]*(-Eqp + Vq) + Iq; */ 549 row[0] = idx + 5; 550 col[0] = idx; 551 col[1] = idx + 1; 552 col[2] = idx + 2; 553 val[0] = -Zdq_inv[3]; 554 val[1] = -Zdq_inv[2]; 555 val[2] = Zdq_inv[2] * dVd_ddelta + Zdq_inv[3] * dVq_ddelta; 556 col[3] = idx + 5; 557 col[4] = net_start + 2 * gbus[i]; 558 col[5] = net_start + 2 * gbus[i] + 1; 559 val[3] = 1; 560 val[4] = Zdq_inv[2] * dVd_dVr + Zdq_inv[3] * dVq_dVr; 561 val[5] = Zdq_inv[2] * dVd_dVi + Zdq_inv[3] * dVq_dVi; 562 PetscCall(MatSetValues(J, 1, row, 6, col, val, INSERT_VALUES)); 563 564 dIGr_ddelta = Id * PetscCosScalar(delta) - Iq * PetscSinScalar(delta); 565 dIGi_ddelta = Id * PetscSinScalar(delta) + Iq * PetscCosScalar(delta); 566 dIGr_dId = PetscSinScalar(delta); 567 dIGr_dIq = PetscCosScalar(delta); 568 dIGi_dId = -PetscCosScalar(delta); 569 dIGi_dIq = PetscSinScalar(delta); 570 571 /* fnet[2*gbus[i]] -= IGi; */ 572 row[0] = net_start + 2 * gbus[i]; 573 col[0] = idx + 2; 574 col[1] = idx + 4; 575 col[2] = idx + 5; 576 val[0] = -dIGi_ddelta; 577 val[1] = -dIGi_dId; 578 val[2] = -dIGi_dIq; 579 PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES)); 580 581 /* fnet[2*gbus[i]+1] -= IGr; */ 582 row[0] = net_start + 2 * gbus[i] + 1; 583 col[0] = idx + 2; 584 col[1] = idx + 4; 585 col[2] = idx + 5; 586 val[0] = -dIGr_ddelta; 587 val[1] = -dIGr_dId; 588 val[2] = -dIGr_dIq; 589 PetscCall(MatSetValues(J, 1, row, 3, col, val, INSERT_VALUES)); 590 591 Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq); 592 593 /* fgen[idx+6] = (KE[i]*Efd + SE - VR)/TE[i]; */ 594 /* SE = k1[i]*PetscExpScalar(k2[i]*Efd); */ 595 dSE_dEfd = k1[i] * k2[i] * PetscExpScalar(k2[i] * Efd); 596 597 row[0] = idx + 6; 598 col[0] = idx + 6; 599 col[1] = idx + 8; 600 val[0] = (KE[i] + dSE_dEfd) / TE[i]; 601 val[1] = -1 / TE[i]; 602 PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES)); 603 604 /* Exciter differential equations */ 605 606 /* fgen[idx+7] = (RF - KF[i]*Efd/TF[i])/TF[i]; */ 607 row[0] = idx + 7; 608 col[0] = idx + 6; 609 col[1] = idx + 7; 610 val[0] = (-KF[i] / TF[i]) / TF[i]; 611 val[1] = 1 / TF[i]; 612 PetscCall(MatSetValues(J, 1, row, 2, col, val, INSERT_VALUES)); 613 614 /* fgen[idx+8] = (VR - KA[i]*RF + KA[i]*KF[i]*Efd/TF[i] - KA[i]*(Vref[i] - Vm))/TA[i]; */ 615 /* Vm = (Vd^2 + Vq^2)^0.5; */ 616 dVm_dVd = Vd / Vm; 617 dVm_dVq = Vq / Vm; 618 dVm_dVr = dVm_dVd * dVd_dVr + dVm_dVq * dVq_dVr; 619 dVm_dVi = dVm_dVd * dVd_dVi + dVm_dVq * dVq_dVi; 620 row[0] = idx + 8; 621 col[0] = idx + 6; 622 col[1] = idx + 7; 623 col[2] = idx + 8; 624 val[0] = (KA[i] * KF[i] / TF[i]) / TA[i]; 625 val[1] = -KA[i] / TA[i]; 626 val[2] = 1 / TA[i]; 627 col[3] = net_start + 2 * gbus[i]; 628 col[4] = net_start + 2 * gbus[i] + 1; 629 val[3] = KA[i] * dVm_dVr / TA[i]; 630 val[4] = KA[i] * dVm_dVi / TA[i]; 631 PetscCall(MatSetValues(J, 1, row, 5, col, val, INSERT_VALUES)); 632 idx = idx + 9; 633 } 634 635 for (i = 0; i < nbus; i++) { 636 PetscCall(MatGetRow(user->Ybus, 2 * i, &ncols, &cols, &yvals)); 637 row[0] = net_start + 2 * i; 638 for (k = 0; k < ncols; k++) { 639 col[k] = net_start + cols[k]; 640 val[k] = yvals[k]; 641 } 642 PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES)); 643 PetscCall(MatRestoreRow(user->Ybus, 2 * i, &ncols, &cols, &yvals)); 644 645 PetscCall(MatGetRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals)); 646 row[0] = net_start + 2 * i + 1; 647 for (k = 0; k < ncols; k++) { 648 col[k] = net_start + cols[k]; 649 val[k] = yvals[k]; 650 } 651 PetscCall(MatSetValues(J, 1, row, ncols, col, val, INSERT_VALUES)); 652 PetscCall(MatRestoreRow(user->Ybus, 2 * i + 1, &ncols, &cols, &yvals)); 653 } 654 655 PetscCall(MatAssemblyBegin(J, MAT_FLUSH_ASSEMBLY)); 656 PetscCall(MatAssemblyEnd(J, MAT_FLUSH_ASSEMBLY)); 657 658 PetscCall(VecGetArray(user->V0, &v0)); 659 for (i = 0; i < nload; i++) { 660 Vr = xnet[2 * lbus[i]]; /* Real part of load bus voltage */ 661 Vi = xnet[2 * lbus[i] + 1]; /* Imaginary part of the load bus voltage */ 662 Vm = PetscSqrtScalar(Vr * Vr + Vi * Vi); 663 Vm2 = Vm * Vm; 664 Vm4 = Vm2 * Vm2; 665 Vm0 = PetscSqrtScalar(v0[2 * lbus[i]] * v0[2 * lbus[i]] + v0[2 * lbus[i] + 1] * v0[2 * lbus[i] + 1]); 666 PD = QD = 0.0; 667 dPD_dVr = dPD_dVi = dQD_dVr = dQD_dVi = 0.0; 668 for (k = 0; k < ld_nsegsp[i]; k++) { 669 PD += ld_alphap[k] * PD0[i] * PetscPowScalar((Vm / Vm0), ld_betap[k]); 670 dPD_dVr += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar((1 / Vm0), ld_betap[k]) * Vr * PetscPowScalar(Vm, (ld_betap[k] - 2)); 671 dPD_dVi += ld_alphap[k] * ld_betap[k] * PD0[i] * PetscPowScalar((1 / Vm0), ld_betap[k]) * Vi * PetscPowScalar(Vm, (ld_betap[k] - 2)); 672 } 673 for (k = 0; k < ld_nsegsq[i]; k++) { 674 QD += ld_alphaq[k] * QD0[i] * PetscPowScalar((Vm / Vm0), ld_betaq[k]); 675 dQD_dVr += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar((1 / Vm0), ld_betaq[k]) * Vr * PetscPowScalar(Vm, (ld_betaq[k] - 2)); 676 dQD_dVi += ld_alphaq[k] * ld_betaq[k] * QD0[i] * PetscPowScalar((1 / Vm0), ld_betaq[k]) * Vi * PetscPowScalar(Vm, (ld_betaq[k] - 2)); 677 } 678 679 /* IDr = (PD*Vr + QD*Vi)/Vm2; */ 680 /* IDi = (-QD*Vr + PD*Vi)/Vm2; */ 681 682 dIDr_dVr = (dPD_dVr * Vr + dQD_dVr * Vi + PD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vr) / Vm4; 683 dIDr_dVi = (dPD_dVi * Vr + dQD_dVi * Vi + QD) / Vm2 - ((PD * Vr + QD * Vi) * 2 * Vi) / Vm4; 684 685 dIDi_dVr = (-dQD_dVr * Vr + dPD_dVr * Vi - QD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vr) / Vm4; 686 dIDi_dVi = (-dQD_dVi * Vr + dPD_dVi * Vi + PD) / Vm2 - ((-QD * Vr + PD * Vi) * 2 * Vi) / Vm4; 687 688 /* fnet[2*lbus[i]] += IDi; */ 689 row[0] = net_start + 2 * lbus[i]; 690 col[0] = net_start + 2 * lbus[i]; 691 col[1] = net_start + 2 * lbus[i] + 1; 692 val[0] = dIDi_dVr; 693 val[1] = dIDi_dVi; 694 PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES)); 695 /* fnet[2*lbus[i]+1] += IDr; */ 696 row[0] = net_start + 2 * lbus[i] + 1; 697 col[0] = net_start + 2 * lbus[i]; 698 col[1] = net_start + 2 * lbus[i] + 1; 699 val[0] = dIDr_dVr; 700 val[1] = dIDr_dVi; 701 PetscCall(MatSetValues(J, 1, row, 2, col, val, ADD_VALUES)); 702 } 703 PetscCall(VecRestoreArray(user->V0, &v0)); 704 705 PetscCall(VecRestoreArray(Xgen, &xgen)); 706 PetscCall(VecRestoreArray(Xnet, &xnet)); 707 708 PetscCall(DMCompositeRestoreLocalVectors(user->dmpgrid, &Xgen, &Xnet)); 709 710 PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY)); 711 PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY)); 712 PetscFunctionReturn(PETSC_SUCCESS); 713 } 714 715 /* 716 J = [I, 0 717 dg_dx, dg_dy] 718 */ 719 PetscErrorCode AlgJacobian(SNES snes, Vec X, Mat A, Mat B, void *ctx) 720 { 721 Userctx *user = (Userctx *)ctx; 722 723 PetscFunctionBegin; 724 PetscCall(ResidualJacobian(snes, X, A, B, ctx)); 725 PetscCall(MatSetOption(A, MAT_KEEP_NONZERO_PATTERN, PETSC_TRUE)); 726 PetscCall(MatZeroRowsIS(A, user->is_diff, 1.0, NULL, NULL)); 727 PetscFunctionReturn(PETSC_SUCCESS); 728 } 729 730 /* 731 J = [a*I-df_dx, -df_dy 732 dg_dx, dg_dy] 733 */ 734 735 PetscErrorCode IJacobian(TS ts, PetscReal t, Vec X, Vec Xdot, PetscReal a, Mat A, Mat B, Userctx *user) 736 { 737 SNES snes; 738 PetscScalar atmp = (PetscScalar)a; 739 PetscInt i, row; 740 741 PetscFunctionBegin; 742 user->t = t; 743 744 PetscCall(TSGetSNES(ts, &snes)); 745 PetscCall(ResidualJacobian(snes, X, A, B, user)); 746 for (i = 0; i < ngen; i++) { 747 row = 9 * i; 748 PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES)); 749 row = 9 * i + 1; 750 PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES)); 751 row = 9 * i + 2; 752 PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES)); 753 row = 9 * i + 3; 754 PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES)); 755 row = 9 * i + 6; 756 PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES)); 757 row = 9 * i + 7; 758 PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES)); 759 row = 9 * i + 8; 760 PetscCall(MatSetValues(A, 1, &row, 1, &row, &atmp, ADD_VALUES)); 761 } 762 PetscCall(MatAssemblyBegin(A, MAT_FINAL_ASSEMBLY)); 763 PetscCall(MatAssemblyEnd(A, MAT_FINAL_ASSEMBLY)); 764 PetscFunctionReturn(PETSC_SUCCESS); 765 } 766 767 int main(int argc, char **argv) 768 { 769 TS ts; 770 SNES snes_alg; 771 PetscMPIInt size; 772 Userctx user; 773 PetscViewer Xview, Ybusview; 774 Vec X; 775 Mat J; 776 PetscInt i; 777 /* sensitivity context */ 778 PetscScalar *y_ptr; 779 Vec lambda[1]; 780 PetscInt *idx2; 781 Vec Xdot; 782 Vec F_alg; 783 PetscInt row_loc, col_loc; 784 PetscScalar val; 785 786 PetscFunctionBeginUser; 787 PetscCall(PetscInitialize(&argc, &argv, "petscoptions", help)); 788 PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size)); 789 PetscCheck(size == 1, PETSC_COMM_WORLD, PETSC_ERR_WRONG_MPI_SIZE, "Only for sequential runs"); 790 791 user.neqs_gen = 9 * ngen; /* # eqs. for generator subsystem */ 792 user.neqs_net = 2 * nbus; /* # eqs. for network subsystem */ 793 user.neqs_pgrid = user.neqs_gen + user.neqs_net; 794 795 /* Create indices for differential and algebraic equations */ 796 PetscCall(PetscMalloc1(7 * ngen, &idx2)); 797 for (i = 0; i < ngen; i++) { 798 idx2[7 * i] = 9 * i; 799 idx2[7 * i + 1] = 9 * i + 1; 800 idx2[7 * i + 2] = 9 * i + 2; 801 idx2[7 * i + 3] = 9 * i + 3; 802 idx2[7 * i + 4] = 9 * i + 6; 803 idx2[7 * i + 5] = 9 * i + 7; 804 idx2[7 * i + 6] = 9 * i + 8; 805 } 806 PetscCall(ISCreateGeneral(PETSC_COMM_WORLD, 7 * ngen, idx2, PETSC_COPY_VALUES, &user.is_diff)); 807 PetscCall(ISComplement(user.is_diff, 0, user.neqs_pgrid, &user.is_alg)); 808 PetscCall(PetscFree(idx2)); 809 810 /* Read initial voltage vector and Ybus */ 811 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "X.bin", FILE_MODE_READ, &Xview)); 812 PetscCall(PetscViewerBinaryOpen(PETSC_COMM_WORLD, "Ybus.bin", FILE_MODE_READ, &Ybusview)); 813 814 PetscCall(VecCreate(PETSC_COMM_WORLD, &user.V0)); 815 PetscCall(VecSetSizes(user.V0, PETSC_DECIDE, user.neqs_net)); 816 PetscCall(VecLoad(user.V0, Xview)); 817 818 PetscCall(MatCreate(PETSC_COMM_WORLD, &user.Ybus)); 819 PetscCall(MatSetSizes(user.Ybus, PETSC_DECIDE, PETSC_DECIDE, user.neqs_net, user.neqs_net)); 820 PetscCall(MatSetType(user.Ybus, MATBAIJ)); 821 /* PetscCall(MatSetBlockSize(user.Ybus,2)); */ 822 PetscCall(MatLoad(user.Ybus, Ybusview)); 823 824 /* Set run time options */ 825 PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", ""); 826 { 827 user.tfaulton = 1.0; 828 user.tfaultoff = 1.2; 829 user.Rfault = 0.0001; 830 user.faultbus = 8; 831 PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL)); 832 PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL)); 833 PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL)); 834 user.t0 = 0.0; 835 user.tmax = 5.0; 836 PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL)); 837 PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL)); 838 } 839 PetscOptionsEnd(); 840 841 PetscCall(PetscViewerDestroy(&Xview)); 842 PetscCall(PetscViewerDestroy(&Ybusview)); 843 844 /* Create DMs for generator and network subsystems */ 845 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_gen, 1, 1, NULL, &user.dmgen)); 846 PetscCall(DMSetOptionsPrefix(user.dmgen, "dmgen_")); 847 PetscCall(DMSetFromOptions(user.dmgen)); 848 PetscCall(DMSetUp(user.dmgen)); 849 PetscCall(DMDACreate1d(PETSC_COMM_WORLD, DM_BOUNDARY_NONE, user.neqs_net, 1, 1, NULL, &user.dmnet)); 850 PetscCall(DMSetOptionsPrefix(user.dmnet, "dmnet_")); 851 PetscCall(DMSetFromOptions(user.dmnet)); 852 PetscCall(DMSetUp(user.dmnet)); 853 /* Create a composite DM packer and add the two DMs */ 854 PetscCall(DMCompositeCreate(PETSC_COMM_WORLD, &user.dmpgrid)); 855 PetscCall(DMSetOptionsPrefix(user.dmpgrid, "pgrid_")); 856 PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmgen)); 857 PetscCall(DMCompositeAddDM(user.dmpgrid, user.dmnet)); 858 859 PetscCall(DMCreateGlobalVector(user.dmpgrid, &X)); 860 861 PetscCall(MatCreate(PETSC_COMM_WORLD, &J)); 862 PetscCall(MatSetSizes(J, PETSC_DECIDE, PETSC_DECIDE, user.neqs_pgrid, user.neqs_pgrid)); 863 PetscCall(MatSetFromOptions(J)); 864 PetscCall(PreallocateJacobian(J, &user)); 865 866 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 867 Create timestepping solver context 868 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 869 PetscCall(TSCreate(PETSC_COMM_WORLD, &ts)); 870 PetscCall(TSSetProblemType(ts, TS_NONLINEAR)); 871 PetscCall(TSSetType(ts, TSCN)); 872 PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)IFunction, &user)); 873 PetscCall(TSSetIJacobian(ts, J, J, (TSIJacobian)IJacobian, &user)); 874 PetscCall(TSSetApplicationContext(ts, &user)); 875 876 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 877 Set initial conditions 878 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 879 PetscCall(SetInitialGuess(X, &user)); 880 /* Just to set up the Jacobian structure */ 881 PetscCall(VecDuplicate(X, &Xdot)); 882 PetscCall(IJacobian(ts, 0.0, X, Xdot, 0.0, J, J, &user)); 883 PetscCall(VecDestroy(&Xdot)); 884 885 /* 886 Save trajectory of solution so that TSAdjointSolve() may be used 887 */ 888 PetscCall(TSSetSaveTrajectory(ts)); 889 890 PetscCall(TSSetMaxTime(ts, user.tmax)); 891 PetscCall(TSSetTimeStep(ts, 0.01)); 892 PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_MATCHSTEP)); 893 PetscCall(TSSetFromOptions(ts)); 894 895 user.alg_flg = PETSC_FALSE; 896 /* Prefault period */ 897 PetscCall(TSSolve(ts, X)); 898 899 /* Create the nonlinear solver for solving the algebraic system */ 900 /* Note that although the algebraic system needs to be solved only for 901 Idq and V, we reuse the entire system including xgen. The xgen 902 variables are held constant by setting their residuals to 0 and 903 putting a 1 on the Jacobian diagonal for xgen rows 904 */ 905 PetscCall(VecDuplicate(X, &F_alg)); 906 PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg)); 907 PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, &user)); 908 PetscCall(MatZeroEntries(J)); 909 PetscCall(SNESSetJacobian(snes_alg, J, J, AlgJacobian, &user)); 910 PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_")); 911 PetscCall(SNESSetFromOptions(snes_alg)); 912 913 /* Apply disturbance - resistive fault at user.faultbus */ 914 /* This is done by adding shunt conductance to the diagonal location 915 in the Ybus matrix */ 916 row_loc = 2 * user.faultbus; 917 col_loc = 2 * user.faultbus + 1; /* Location for G */ 918 val = 1 / user.Rfault; 919 PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES)); 920 row_loc = 2 * user.faultbus + 1; 921 col_loc = 2 * user.faultbus; /* Location for G */ 922 val = 1 / user.Rfault; 923 PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES)); 924 925 PetscCall(MatAssemblyBegin(user.Ybus, MAT_FINAL_ASSEMBLY)); 926 PetscCall(MatAssemblyEnd(user.Ybus, MAT_FINAL_ASSEMBLY)); 927 928 user.alg_flg = PETSC_TRUE; 929 /* Solve the algebraic equations */ 930 PetscCall(SNESSolve(snes_alg, NULL, X)); 931 932 /* Disturbance period */ 933 user.alg_flg = PETSC_FALSE; 934 PetscCall(TSSetTime(ts, user.tfaulton)); 935 PetscCall(TSSetMaxTime(ts, user.tfaultoff)); 936 PetscCall(TSSolve(ts, X)); 937 938 /* Remove the fault */ 939 row_loc = 2 * user.faultbus; 940 col_loc = 2 * user.faultbus + 1; 941 val = -1 / user.Rfault; 942 PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES)); 943 row_loc = 2 * user.faultbus + 1; 944 col_loc = 2 * user.faultbus; 945 val = -1 / user.Rfault; 946 PetscCall(MatSetValues(user.Ybus, 1, &row_loc, 1, &col_loc, &val, ADD_VALUES)); 947 948 PetscCall(MatAssemblyBegin(user.Ybus, MAT_FINAL_ASSEMBLY)); 949 PetscCall(MatAssemblyEnd(user.Ybus, MAT_FINAL_ASSEMBLY)); 950 951 PetscCall(MatZeroEntries(J)); 952 953 user.alg_flg = PETSC_TRUE; 954 955 /* Solve the algebraic equations */ 956 PetscCall(SNESSolve(snes_alg, NULL, X)); 957 958 /* Post-disturbance period */ 959 user.alg_flg = PETSC_TRUE; 960 PetscCall(TSSetTime(ts, user.tfaultoff)); 961 PetscCall(TSSetMaxTime(ts, user.tmax)); 962 PetscCall(TSSolve(ts, X)); 963 964 /* - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - 965 Adjoint model starts here 966 - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - - */ 967 PetscCall(TSSetPostStep(ts, NULL)); 968 PetscCall(MatCreateVecs(J, &lambda[0], NULL)); 969 /* Set initial conditions for the adjoint integration */ 970 PetscCall(VecZeroEntries(lambda[0])); 971 PetscCall(VecGetArray(lambda[0], &y_ptr)); 972 y_ptr[0] = 1.0; 973 PetscCall(VecRestoreArray(lambda[0], &y_ptr)); 974 PetscCall(TSSetCostGradients(ts, 1, lambda, NULL)); 975 976 PetscCall(TSAdjointSolve(ts)); 977 978 PetscCall(PetscPrintf(PETSC_COMM_WORLD, "\n sensitivity wrt initial conditions: \n")); 979 PetscCall(VecView(lambda[0], PETSC_VIEWER_STDOUT_WORLD)); 980 PetscCall(VecDestroy(&lambda[0])); 981 982 PetscCall(SNESDestroy(&snes_alg)); 983 PetscCall(VecDestroy(&F_alg)); 984 PetscCall(MatDestroy(&J)); 985 PetscCall(MatDestroy(&user.Ybus)); 986 PetscCall(VecDestroy(&X)); 987 PetscCall(VecDestroy(&user.V0)); 988 PetscCall(DMDestroy(&user.dmgen)); 989 PetscCall(DMDestroy(&user.dmnet)); 990 PetscCall(DMDestroy(&user.dmpgrid)); 991 PetscCall(ISDestroy(&user.is_diff)); 992 PetscCall(ISDestroy(&user.is_alg)); 993 PetscCall(TSDestroy(&ts)); 994 PetscCall(PetscFinalize()); 995 return 0; 996 } 997 998 /*TEST 999 1000 build: 1001 requires: double !complex !defined(PETSC_USE_64BIT_INDICES) 1002 1003 test: 1004 args: -viewer_binary_skip_info 1005 localrunfiles: petscoptions X.bin Ybus.bin 1006 1007 test: 1008 suffix: 2 1009 args: -viewer_binary_skip_info -ts_type beuler 1010 localrunfiles: petscoptions X.bin Ybus.bin 1011 1012 TEST*/ 1013