xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9busdmnetwork.c (revision 58d68138c660dfb4e9f5b03334792cd4f2ffd7cc)
1 
2 static char help[] = "This example uses the same problem set up of ex9busdmnetwork.c. \n\
3 It demonstrates setting and accessing of variables for individual components, instead of \n\
4 the network vertices (as used in ex9busdmnetwork.c). This is especially useful where vertices \n\
5 /edges have multiple-components associated with them and one or more components has physics \n\
6 associated with it. \n\
7 Input parameters include:\n\
8   -nc : number of copies of the base case\n\n";
9 
10 /*
11    This example was modified from ex9busdmnetwork.c.
12 */
13 
14 #include <petscts.h>
15 #include <petscdmnetwork.h>
16 
17 #define FREQ    60
18 #define W_S     (2 * PETSC_PI * FREQ)
19 #define NGEN    3 /* No. of generators in the 9 bus system */
20 #define NLOAD   3 /* No. of loads in the 9 bus system */
21 #define NBUS    9 /* No. of buses in the 9 bus system */
22 #define NBRANCH 9 /* No. of branches in the 9 bus system */
23 
24 typedef struct {
25   PetscInt    id;    /* Bus Number or extended bus name*/
26   PetscScalar mbase; /* MVA base of the machine */
27   PetscScalar PG;    /* Generator active power output */
28   PetscScalar QG;    /* Generator reactive power output */
29 
30   /* Generator constants */
31   PetscScalar H;    /* Inertia constant */
32   PetscScalar Rs;   /* Stator Resistance */
33   PetscScalar Xd;   /* d-axis reactance */
34   PetscScalar Xdp;  /* d-axis transient reactance */
35   PetscScalar Xq;   /* q-axis reactance Xq(1) set to 0.4360, value given in text 0.0969 */
36   PetscScalar Xqp;  /* q-axis transient reactance */
37   PetscScalar Td0p; /* d-axis open circuit time constant */
38   PetscScalar Tq0p; /* q-axis open circuit time constant */
39   PetscScalar M;    /* M = 2*H/W_S */
40   PetscScalar D;    /* D = 0.1*M */
41   PetscScalar TM;   /* Mechanical Torque */
42 } Gen;
43 
44 typedef struct {
45   /* Exciter system constants */
46   PetscScalar KA;     /* Voltage regulator gain constant */
47   PetscScalar TA;     /* Voltage regulator time constant */
48   PetscScalar KE;     /* Exciter gain constant */
49   PetscScalar TE;     /* Exciter time constant */
50   PetscScalar KF;     /* Feedback stabilizer gain constant */
51   PetscScalar TF;     /* Feedback stabilizer time constant */
52   PetscScalar k1, k2; /* calculating the saturation function SE = k1*exp(k2*Efd) */
53   PetscScalar Vref;   /* Voltage regulator voltage setpoint */
54 } Exc;
55 
56 typedef struct {
57   PetscInt    id;      /* node id */
58   PetscInt    nofgen;  /* Number of generators at the bus*/
59   PetscInt    nofload; /*  Number of load at the bus*/
60   PetscScalar yff[2];  /* yff[0]= imaginary part of admittance, yff[1]=real part of admittance*/
61   PetscScalar vr;      /* Real component of bus voltage */
62   PetscScalar vi;      /* Imaginary component of bus voltage */
63 } Bus;
64 
65 /* Load constants
66   We use a composite load model that describes the load and reactive powers at each time instant as follows
67   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
68   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
69   where
70     id                  - index of the load
71     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
72     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
73     P_D0                - Real power load
74     Q_D0                - Reactive power load
75     Vm(t)              - Voltage magnitude at time t
76     Vm0                - Voltage magnitude at t = 0
77     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
78 
79     Note: All loads have the same characteristic currently.
80   */
81 typedef struct {
82   PetscInt    id; /* bus id */
83   PetscInt    ld_nsegsp, ld_nsegsq;
84   PetscScalar PD0, QD0;
85   PetscScalar ld_alphap[3]; /* ld_alphap=[1,0,0], an array, not a value, so use *ld_alphap; */
86   PetscScalar ld_betap[3], ld_alphaq[3], ld_betaq[3];
87 } Load;
88 
89 typedef struct {
90   PetscInt    id;     /* node id */
91   PetscScalar yft[2]; /* yft[0]= imaginary part of admittance, yft[1]=real part of admittance*/
92 } Branch;
93 
94 typedef struct {
95   PetscReal    tfaulton, tfaultoff; /* Fault on and off times */
96   PetscReal    t;
97   PetscReal    t0, tmax; /* initial time and final time */
98   PetscInt     faultbus; /* Fault bus */
99   PetscScalar  Rfault;   /* Fault resistance (pu) */
100   PetscScalar *ybusfault;
101   PetscBool    alg_flg;
102 } Userctx;
103 
104 /* Used to read data into the DMNetwork components */
105 PetscErrorCode read_data(PetscInt nc, Gen **pgen, Exc **pexc, Load **pload, Bus **pbus, Branch **pbranch, PetscInt **pedgelist) {
106   PetscInt           i, j, row[1], col[2];
107   PetscInt          *edgelist;
108   PetscInt           nofgen[9]  = {1, 1, 1, 0, 0, 0, 0, 0, 0}; /* Buses at which generators are incident */
109   PetscInt           nofload[9] = {0, 0, 0, 0, 1, 1, 0, 1, 0}; /* Buses at which loads are incident */
110   const PetscScalar *varr;
111   PetscScalar        M[3], D[3];
112   Bus               *bus;
113   Branch            *branch;
114   Gen               *gen;
115   Exc               *exc;
116   Load              *load;
117   Mat                Ybus;
118   Vec                V0;
119 
120   /*10 parameters*/
121   /* Generator real and reactive powers (found via loadflow) */
122   static const PetscScalar PG[3] = {0.716786142395021, 1.630000000000000, 0.850000000000000};
123   static const PetscScalar QG[3] = {0.270702180178785, 0.066120127797275, -0.108402221791588};
124 
125   /* Generator constants */
126   static const PetscScalar H[3]    = {23.64, 6.4, 3.01};       /* Inertia constant */
127   static const PetscScalar Rs[3]   = {0.0, 0.0, 0.0};          /* Stator Resistance */
128   static const PetscScalar Xd[3]   = {0.146, 0.8958, 1.3125};  /* d-axis reactance */
129   static const PetscScalar Xdp[3]  = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
130   static 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 */
131   static const PetscScalar Xqp[3]  = {0.0969, 0.1969, 0.25};   /* q-axis transient reactance */
132   static const PetscScalar Td0p[3] = {8.96, 6.0, 5.89};        /* d-axis open circuit time constant */
133   static const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6};       /* q-axis open circuit time constant */
134 
135   /* Exciter system constants (8 parameters)*/
136   static const PetscScalar KA[3] = {20.0, 20.0, 20.0};    /* Voltage regulartor gain constant */
137   static const PetscScalar TA[3] = {0.2, 0.2, 0.2};       /* Voltage regulator time constant */
138   static const PetscScalar KE[3] = {1.0, 1.0, 1.0};       /* Exciter gain constant */
139   static const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */
140   static const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
141   static const PetscScalar TF[3] = {0.35, 0.35, 0.35};    /* Feedback stabilizer time constant */
142   static const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039};
143   static const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
144 
145   /* Load constants */
146   static const PetscScalar PD0[3]       = {1.25, 0.9, 1.0};
147   static const PetscScalar QD0[3]       = {0.5, 0.3, 0.35};
148   static const PetscScalar ld_alphaq[3] = {1, 0, 0};
149   static const PetscScalar ld_betaq[3]  = {2, 1, 0};
150   static const PetscScalar ld_betap[3]  = {2, 1, 0};
151   static const PetscScalar ld_alphap[3] = {1, 0, 0};
152   PetscInt                 ld_nsegsp[3] = {3, 3, 3};
153   PetscInt                 ld_nsegsq[3] = {3, 3, 3};
154   PetscViewer              Xview, Ybusview;
155   PetscInt                 neqs_net, m, n;
156 
157   PetscFunctionBeginUser;
158   /* Read V0 and Ybus from files */
159   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "X.bin", FILE_MODE_READ, &Xview));
160   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "Ybus.bin", FILE_MODE_READ, &Ybusview));
161   PetscCall(VecCreate(PETSC_COMM_SELF, &V0));
162   PetscCall(VecLoad(V0, Xview));
163 
164   PetscCall(MatCreate(PETSC_COMM_SELF, &Ybus));
165   PetscCall(MatSetType(Ybus, MATBAIJ));
166   PetscCall(MatLoad(Ybus, Ybusview));
167 
168   /* Destroy unnecessary stuff */
169   PetscCall(PetscViewerDestroy(&Xview));
170   PetscCall(PetscViewerDestroy(&Ybusview));
171 
172   PetscCall(MatGetLocalSize(Ybus, &m, &n));
173   neqs_net = 2 * NBUS; /* # eqs. for network subsystem   */
174   PetscCheck(m == neqs_net && n == neqs_net, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matrix Ybus is in wrong sizes");
175 
176   M[0] = 2 * H[0] / W_S;
177   M[1] = 2 * H[1] / W_S;
178   M[2] = 2 * H[2] / W_S;
179   D[0] = 0.1 * M[0];
180   D[1] = 0.1 * M[1];
181   D[2] = 0.1 * M[2];
182 
183   /* Alocate memory for bus, generators, exciter, loads and branches */
184   PetscCall(PetscCalloc5(NBUS * nc, &bus, NGEN * nc, &gen, NLOAD * nc, &load, NBRANCH * nc + (nc - 1), &branch, NGEN * nc, &exc));
185 
186   PetscCall(VecGetArrayRead(V0, &varr));
187 
188   /* read bus data */
189   for (i = 0; i < nc; i++) {
190     for (j = 0; j < NBUS; j++) {
191       bus[i * 9 + j].id      = i * 9 + j;
192       bus[i * 9 + j].nofgen  = nofgen[j];
193       bus[i * 9 + j].nofload = nofload[j];
194       bus[i * 9 + j].vr      = varr[2 * j];
195       bus[i * 9 + j].vi      = varr[2 * j + 1];
196       row[0]                 = 2 * j;
197       col[0]                 = 2 * j;
198       col[1]                 = 2 * j + 1;
199       /* real and imaginary part of admittance from Ybus into yff */
200       PetscCall(MatGetValues(Ybus, 1, row, 2, col, bus[i * 9 + j].yff));
201     }
202   }
203 
204   /* read generator data */
205   for (i = 0; i < nc; i++) {
206     for (j = 0; j < NGEN; j++) {
207       gen[i * 3 + j].id   = i * 3 + j;
208       gen[i * 3 + j].PG   = PG[j];
209       gen[i * 3 + j].QG   = QG[j];
210       gen[i * 3 + j].H    = H[j];
211       gen[i * 3 + j].Rs   = Rs[j];
212       gen[i * 3 + j].Xd   = Xd[j];
213       gen[i * 3 + j].Xdp  = Xdp[j];
214       gen[i * 3 + j].Xq   = Xq[j];
215       gen[i * 3 + j].Xqp  = Xqp[j];
216       gen[i * 3 + j].Td0p = Td0p[j];
217       gen[i * 3 + j].Tq0p = Tq0p[j];
218       gen[i * 3 + j].M    = M[j];
219       gen[i * 3 + j].D    = D[j];
220     }
221   }
222 
223   for (i = 0; i < nc; i++) {
224     for (j = 0; j < NGEN; j++) {
225       /* exciter system */
226       exc[i * 3 + j].KA = KA[j];
227       exc[i * 3 + j].TA = TA[j];
228       exc[i * 3 + j].KE = KE[j];
229       exc[i * 3 + j].TE = TE[j];
230       exc[i * 3 + j].KF = KF[j];
231       exc[i * 3 + j].TF = TF[j];
232       exc[i * 3 + j].k1 = k1[j];
233       exc[i * 3 + j].k2 = k2[j];
234     }
235   }
236 
237   /* read load data */
238   for (i = 0; i < nc; i++) {
239     for (j = 0; j < NLOAD; j++) {
240       load[i * 3 + j].id        = i * 3 + j;
241       load[i * 3 + j].PD0       = PD0[j];
242       load[i * 3 + j].QD0       = QD0[j];
243       load[i * 3 + j].ld_nsegsp = ld_nsegsp[j];
244 
245       load[i * 3 + j].ld_alphap[0] = ld_alphap[0];
246       load[i * 3 + j].ld_alphap[1] = ld_alphap[1];
247       load[i * 3 + j].ld_alphap[2] = ld_alphap[2];
248 
249       load[i * 3 + j].ld_alphaq[0] = ld_alphaq[0];
250       load[i * 3 + j].ld_alphaq[1] = ld_alphaq[1];
251       load[i * 3 + j].ld_alphaq[2] = ld_alphaq[2];
252 
253       load[i * 3 + j].ld_betap[0] = ld_betap[0];
254       load[i * 3 + j].ld_betap[1] = ld_betap[1];
255       load[i * 3 + j].ld_betap[2] = ld_betap[2];
256       load[i * 3 + j].ld_nsegsq   = ld_nsegsq[j];
257 
258       load[i * 3 + j].ld_betaq[0] = ld_betaq[0];
259       load[i * 3 + j].ld_betaq[1] = ld_betaq[1];
260       load[i * 3 + j].ld_betaq[2] = ld_betaq[2];
261     }
262   }
263   PetscCall(PetscCalloc1(2 * NBRANCH * nc + 2 * (nc - 1), &edgelist));
264 
265   /* read edgelist */
266   for (i = 0; i < nc; i++) {
267     for (j = 0; j < NBRANCH; j++) {
268       switch (j) {
269       case 0:
270         edgelist[i * 18 + 2 * j]     = 0 + 9 * i;
271         edgelist[i * 18 + 2 * j + 1] = 3 + 9 * i;
272         break;
273       case 1:
274         edgelist[i * 18 + 2 * j]     = 1 + 9 * i;
275         edgelist[i * 18 + 2 * j + 1] = 6 + 9 * i;
276         break;
277       case 2:
278         edgelist[i * 18 + 2 * j]     = 2 + 9 * i;
279         edgelist[i * 18 + 2 * j + 1] = 8 + 9 * i;
280         break;
281       case 3:
282         edgelist[i * 18 + 2 * j]     = 3 + 9 * i;
283         edgelist[i * 18 + 2 * j + 1] = 4 + 9 * i;
284         break;
285       case 4:
286         edgelist[i * 18 + 2 * j]     = 3 + 9 * i;
287         edgelist[i * 18 + 2 * j + 1] = 5 + 9 * i;
288         break;
289       case 5:
290         edgelist[i * 18 + 2 * j]     = 4 + 9 * i;
291         edgelist[i * 18 + 2 * j + 1] = 6 + 9 * i;
292         break;
293       case 6:
294         edgelist[i * 18 + 2 * j]     = 5 + 9 * i;
295         edgelist[i * 18 + 2 * j + 1] = 8 + 9 * i;
296         break;
297       case 7:
298         edgelist[i * 18 + 2 * j]     = 6 + 9 * i;
299         edgelist[i * 18 + 2 * j + 1] = 7 + 9 * i;
300         break;
301       case 8:
302         edgelist[i * 18 + 2 * j]     = 7 + 9 * i;
303         edgelist[i * 18 + 2 * j + 1] = 8 + 9 * i;
304         break;
305       default: break;
306       }
307     }
308   }
309 
310   /* for connecting last bus of previous network(9*i-1) to first bus of next network(9*i), the branch admittance=-0.0301407+j17.3611 */
311   for (i = 1; i < nc; i++) {
312     edgelist[18 * nc + 2 * (i - 1)]     = 8 + (i - 1) * 9;
313     edgelist[18 * nc + 2 * (i - 1) + 1] = 9 * i;
314 
315     /* adding admittances to the off-diagonal elements */
316     branch[9 * nc + (i - 1)].id     = 9 * nc + (i - 1);
317     branch[9 * nc + (i - 1)].yft[0] = 17.3611;
318     branch[9 * nc + (i - 1)].yft[1] = -0.0301407;
319 
320     /* subtracting admittances from the diagonal elements */
321     bus[9 * i - 1].yff[0] -= 17.3611;
322     bus[9 * i - 1].yff[1] -= -0.0301407;
323 
324     bus[9 * i].yff[0] -= 17.3611;
325     bus[9 * i].yff[1] -= -0.0301407;
326   }
327 
328   /* read branch data */
329   for (i = 0; i < nc; i++) {
330     for (j = 0; j < NBRANCH; j++) {
331       branch[i * 9 + j].id = i * 9 + j;
332 
333       row[0] = edgelist[2 * j] * 2;
334       col[0] = edgelist[2 * j + 1] * 2;
335       col[1] = edgelist[2 * j + 1] * 2 + 1;
336       PetscCall(MatGetValues(Ybus, 1, row, 2, col, branch[i * 9 + j].yft)); /*imaginary part of admittance*/
337     }
338   }
339 
340   *pgen      = gen;
341   *pexc      = exc;
342   *pload     = load;
343   *pbus      = bus;
344   *pbranch   = branch;
345   *pedgelist = edgelist;
346 
347   PetscCall(VecRestoreArrayRead(V0, &varr));
348 
349   /* Destroy unnecessary stuff */
350   PetscCall(MatDestroy(&Ybus));
351   PetscCall(VecDestroy(&V0));
352   PetscFunctionReturn(0);
353 }
354 
355 PetscErrorCode SetInitialGuess(DM networkdm, Vec X) {
356   Bus         *bus;
357   Gen         *gen;
358   Exc         *exc;
359   PetscInt     v, vStart, vEnd, offset;
360   PetscInt     key, numComps, j;
361   PetscBool    ghostvtex;
362   Vec          localX;
363   PetscScalar *xarr;
364   PetscScalar  Vr = 0, Vi = 0, Vm = 0, Vm2; /* Terminal voltage variables */
365   PetscScalar  IGr, IGi;                    /* Generator real and imaginary current */
366   PetscScalar  Eqp, Edp, delta;             /* Generator variables */
367   PetscScalar  Efd = 0, RF, VR;             /* Exciter variables */
368   PetscScalar  Vd, Vq;                      /* Generator dq axis voltages */
369   PetscScalar  Id, Iq;                      /* Generator dq axis currents */
370   PetscScalar  theta;                       /* Generator phase angle */
371   PetscScalar  SE;
372   void        *component;
373 
374   PetscFunctionBegin;
375   PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
376   PetscCall(DMGetLocalVector(networkdm, &localX));
377 
378   PetscCall(VecSet(X, 0.0));
379   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
380   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
381 
382   PetscCall(VecGetArray(localX, &xarr));
383 
384   for (v = vStart; v < vEnd; v++) {
385     PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
386     if (ghostvtex) continue;
387 
388     PetscCall(DMNetworkGetNumComponents(networkdm, v, &numComps));
389     for (j = 0; j < numComps; j++) {
390       PetscCall(DMNetworkGetComponent(networkdm, v, j, &key, &component, NULL));
391       if (key == 1) {
392         bus = (Bus *)(component);
393 
394         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offset));
395         xarr[offset]     = bus->vr;
396         xarr[offset + 1] = bus->vi;
397 
398         Vr = bus->vr;
399         Vi = bus->vi;
400       } else if (key == 2) {
401         gen = (Gen *)(component);
402         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offset));
403         Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
404         Vm2 = Vm * Vm;
405         /* Real part of gen current */
406         IGr = (Vr * gen->PG + Vi * gen->QG) / Vm2;
407         /* Imaginary part of gen current */
408         IGi = (Vi * gen->PG - Vr * gen->QG) / Vm2;
409 
410         /* Machine angle */
411         delta = atan2(Vi + gen->Xq * IGr, Vr - gen->Xq * IGi);
412         theta = PETSC_PI / 2.0 - delta;
413 
414         /* d-axis stator current */
415         Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta);
416 
417         /* q-axis stator current */
418         Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta);
419 
420         Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
421         Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
422 
423         /* d-axis transient EMF */
424         Edp = Vd + gen->Rs * Id - gen->Xqp * Iq;
425 
426         /* q-axis transient EMF */
427         Eqp = Vq + gen->Rs * Iq + gen->Xdp * Id;
428 
429         gen->TM = gen->PG;
430 
431         xarr[offset]     = Eqp;
432         xarr[offset + 1] = Edp;
433         xarr[offset + 2] = delta;
434         xarr[offset + 3] = W_S;
435         xarr[offset + 4] = Id;
436         xarr[offset + 5] = Iq;
437 
438         Efd = Eqp + (gen->Xd - gen->Xdp) * Id;
439 
440       } else if (key == 3) {
441         exc = (Exc *)(component);
442         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offset));
443 
444         SE = exc->k1 * PetscExpScalar(exc->k2 * Efd);
445         VR = exc->KE * Efd + SE;
446         RF = exc->KF * Efd / exc->TF;
447 
448         xarr[offset]     = Efd;
449         xarr[offset + 1] = RF;
450         xarr[offset + 2] = VR;
451 
452         exc->Vref = Vm + (VR / exc->KA);
453       }
454     }
455   }
456   PetscCall(VecRestoreArray(localX, &xarr));
457   PetscCall(DMLocalToGlobalBegin(networkdm, localX, ADD_VALUES, X));
458   PetscCall(DMLocalToGlobalEnd(networkdm, localX, ADD_VALUES, X));
459   PetscCall(DMRestoreLocalVector(networkdm, &localX));
460   PetscFunctionReturn(0);
461 }
462 
463 /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
464 PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi) {
465   PetscFunctionBegin;
466   *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
467   *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
468   PetscFunctionReturn(0);
469 }
470 
471 /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
472 PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq) {
473   PetscFunctionBegin;
474   *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
475   *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
476   PetscFunctionReturn(0);
477 }
478 
479 /* Computes F(t,U,U_t) where F() = 0 is the DAE to be solved. */
480 PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user) {
481   DM                 networkdm;
482   Vec                localX, localXdot, localF;
483   PetscInt           vfrom, vto, offsetfrom, offsetto;
484   PetscInt           v, vStart, vEnd, e;
485   PetscScalar       *farr;
486   PetscScalar        Vd = 0, Vq = 0, SE;
487   const PetscScalar *xarr, *xdotarr;
488   void              *component;
489   PetscScalar        Vr = 0, Vi = 0;
490 
491   PetscFunctionBegin;
492   PetscCall(VecSet(F, 0.0));
493 
494   PetscCall(TSGetDM(ts, &networkdm));
495   PetscCall(DMGetLocalVector(networkdm, &localF));
496   PetscCall(DMGetLocalVector(networkdm, &localX));
497   PetscCall(DMGetLocalVector(networkdm, &localXdot));
498   PetscCall(VecSet(localF, 0.0));
499 
500   /* update ghost values of localX and localXdot */
501   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
502   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
503 
504   PetscCall(DMGlobalToLocalBegin(networkdm, Xdot, INSERT_VALUES, localXdot));
505   PetscCall(DMGlobalToLocalEnd(networkdm, Xdot, INSERT_VALUES, localXdot));
506 
507   PetscCall(VecGetArrayRead(localX, &xarr));
508   PetscCall(VecGetArrayRead(localXdot, &xdotarr));
509   PetscCall(VecGetArray(localF, &farr));
510 
511   PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
512 
513   for (v = vStart; v < vEnd; v++) {
514     PetscInt    i, j, offsetbus, offsetgen, offsetexc, key;
515     Bus        *bus;
516     Gen        *gen;
517     Exc        *exc;
518     Load       *load;
519     PetscBool   ghostvtex;
520     PetscInt    numComps;
521     PetscScalar Yffr, Yffi; /* Real and imaginary fault admittances */
522     PetscScalar Vm, Vm2, Vm0;
523     PetscScalar Vr0 = 0, Vi0 = 0;
524     PetscScalar PD, QD;
525 
526     PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
527     PetscCall(DMNetworkGetNumComponents(networkdm, v, &numComps));
528 
529     for (j = 0; j < numComps; j++) {
530       PetscCall(DMNetworkGetComponent(networkdm, v, j, &key, &component, NULL));
531       if (key == 1) {
532         PetscInt        nconnedges;
533         const PetscInt *connedges;
534 
535         bus = (Bus *)(component);
536         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetbus));
537         if (!ghostvtex) {
538           Vr = xarr[offsetbus];
539           Vi = xarr[offsetbus + 1];
540 
541           Yffr = bus->yff[1];
542           Yffi = bus->yff[0];
543 
544           if (user->alg_flg) {
545             Yffr += user->ybusfault[bus->id * 2 + 1];
546             Yffi += user->ybusfault[bus->id * 2];
547           }
548           Vr0 = bus->vr;
549           Vi0 = bus->vi;
550 
551           /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
552            The generator current injection, IG, and load current injection, ID are added later
553            */
554           farr[offsetbus] += Yffi * Vr + Yffr * Vi;     /* imaginary current due to diagonal elements */
555           farr[offsetbus + 1] += Yffr * Vr - Yffi * Vi; /* real current due to diagonal elements */
556         }
557 
558         PetscCall(DMNetworkGetSupportingEdges(networkdm, v, &nconnedges, &connedges));
559 
560         for (i = 0; i < nconnedges; i++) {
561           Branch         *branch;
562           PetscInt        keye;
563           PetscScalar     Yfti, Yftr, Vfr, Vfi, Vtr, Vti;
564           const PetscInt *cone;
565 
566           e = connedges[i];
567           PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, (void **)&branch, NULL));
568 
569           Yfti = branch->yft[0];
570           Yftr = branch->yft[1];
571 
572           PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone));
573 
574           vfrom = cone[0];
575           vto   = cone[1];
576 
577           PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, 0, &offsetfrom));
578           PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, 0, &offsetto));
579 
580           /* From bus and to bus real and imaginary voltages */
581           Vfr = xarr[offsetfrom];
582           Vfi = xarr[offsetfrom + 1];
583           Vtr = xarr[offsetto];
584           Vti = xarr[offsetto + 1];
585 
586           if (vfrom == v) {
587             farr[offsetfrom] += Yftr * Vti + Yfti * Vtr;
588             farr[offsetfrom + 1] += Yftr * Vtr - Yfti * Vti;
589           } else {
590             farr[offsetto] += Yftr * Vfi + Yfti * Vfr;
591             farr[offsetto + 1] += Yftr * Vfr - Yfti * Vfi;
592           }
593         }
594       } else if (key == 2) {
595         if (!ghostvtex) {
596           PetscScalar Eqp, Edp, delta, w; /* Generator variables */
597           PetscScalar Efd;                /* Exciter field voltage */
598           PetscScalar Id, Iq;             /* Generator dq axis currents */
599           PetscScalar IGr, IGi, Zdq_inv[4], det;
600           PetscScalar Xd, Xdp, Td0p, Xq, Xqp, Tq0p, TM, D, M, Rs; /* Generator parameters */
601 
602           gen = (Gen *)(component);
603           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetgen));
604 
605           /* Generator state variables */
606           Eqp   = xarr[offsetgen];
607           Edp   = xarr[offsetgen + 1];
608           delta = xarr[offsetgen + 2];
609           w     = xarr[offsetgen + 3];
610           Id    = xarr[offsetgen + 4];
611           Iq    = xarr[offsetgen + 5];
612 
613           /* Generator parameters */
614           Xd   = gen->Xd;
615           Xdp  = gen->Xdp;
616           Td0p = gen->Td0p;
617           Xq   = gen->Xq;
618           Xqp  = gen->Xqp;
619           Tq0p = gen->Tq0p;
620           TM   = gen->TM;
621           D    = gen->D;
622           M    = gen->M;
623           Rs   = gen->Rs;
624 
625           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, 2, &offsetexc));
626           Efd = xarr[offsetexc];
627 
628           /* Generator differential equations */
629           farr[offsetgen]     = (Eqp + (Xd - Xdp) * Id - Efd) / Td0p + xdotarr[offsetgen];
630           farr[offsetgen + 1] = (Edp - (Xq - Xqp) * Iq) / Tq0p + xdotarr[offsetgen + 1];
631           farr[offsetgen + 2] = -w + W_S + xdotarr[offsetgen + 2];
632           farr[offsetgen + 3] = (-TM + Edp * Id + Eqp * Iq + (Xqp - Xdp) * Id * Iq + D * (w - W_S)) / M + xdotarr[offsetgen + 3];
633 
634           PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
635 
636           /* Algebraic equations for stator currents */
637           det = Rs * Rs + Xdp * Xqp;
638 
639           Zdq_inv[0] = Rs / det;
640           Zdq_inv[1] = Xqp / det;
641           Zdq_inv[2] = -Xdp / det;
642           Zdq_inv[3] = Rs / det;
643 
644           farr[offsetgen + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
645           farr[offsetgen + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
646 
647           PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
648 
649           /* Add generator current injection to network */
650           farr[offsetbus] -= IGi;
651           farr[offsetbus + 1] -= IGr;
652         }
653       } else if (key == 3) {
654         if (!ghostvtex) {
655           PetscScalar k1, k2, KE, TE, TF, KA, KF, Vref, TA; /* Generator parameters */
656           PetscScalar Efd, RF, VR;                          /* Exciter variables */
657 
658           exc = (Exc *)(component);
659           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetexc));
660 
661           Efd = xarr[offsetexc];
662           RF  = xarr[offsetexc + 1];
663           VR  = xarr[offsetexc + 2];
664 
665           k1   = exc->k1;
666           k2   = exc->k2;
667           KE   = exc->KE;
668           TE   = exc->TE;
669           TF   = exc->TF;
670           KA   = exc->KA;
671           KF   = exc->KF;
672           Vref = exc->Vref;
673           TA   = exc->TA;
674 
675           Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
676           SE = k1 * PetscExpScalar(k2 * Efd);
677 
678           /* Exciter differential equations */
679           farr[offsetexc]     = (KE * Efd + SE - VR) / TE + xdotarr[offsetexc];
680           farr[offsetexc + 1] = (RF - KF * Efd / TF) / TF + xdotarr[offsetexc + 1];
681           farr[offsetexc + 2] = (VR - KA * RF + KA * KF * Efd / TF - KA * (Vref - Vm)) / TA + xdotarr[offsetexc + 2];
682         }
683       } else if (key == 4) {
684         if (!ghostvtex) {
685           PetscInt     k;
686           PetscInt     ld_nsegsp;
687           PetscInt     ld_nsegsq;
688           PetscScalar *ld_alphap;
689           PetscScalar *ld_betap, *ld_alphaq, *ld_betaq, PD0, QD0, IDr, IDi;
690 
691           load = (Load *)(component);
692 
693           /* Load Parameters */
694           ld_nsegsp = load->ld_nsegsp;
695           ld_alphap = load->ld_alphap;
696           ld_betap  = load->ld_betap;
697           ld_nsegsq = load->ld_nsegsq;
698           ld_alphaq = load->ld_alphaq;
699           ld_betaq  = load->ld_betaq;
700           PD0       = load->PD0;
701           QD0       = load->QD0;
702 
703           Vr  = xarr[offsetbus];     /* Real part of generator terminal voltage */
704           Vi  = xarr[offsetbus + 1]; /* Imaginary part of the generator terminal voltage */
705           Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
706           Vm2 = Vm * Vm;
707           Vm0 = PetscSqrtScalar(Vr0 * Vr0 + Vi0 * Vi0);
708           PD = QD = 0.0;
709           for (k = 0; k < ld_nsegsp; k++) PD += ld_alphap[k] * PD0 * PetscPowScalar((Vm / Vm0), ld_betap[k]);
710           for (k = 0; k < ld_nsegsq; k++) QD += ld_alphaq[k] * QD0 * PetscPowScalar((Vm / Vm0), ld_betaq[k]);
711 
712           /* Load currents */
713           IDr = (PD * Vr + QD * Vi) / Vm2;
714           IDi = (-QD * Vr + PD * Vi) / Vm2;
715 
716           /* Load current contribution to the network */
717           farr[offsetbus] += IDi;
718           farr[offsetbus + 1] += IDr;
719         }
720       }
721     }
722   }
723 
724   PetscCall(VecRestoreArrayRead(localX, &xarr));
725   PetscCall(VecRestoreArrayRead(localXdot, &xdotarr));
726   PetscCall(VecRestoreArray(localF, &farr));
727   PetscCall(DMRestoreLocalVector(networkdm, &localX));
728   PetscCall(DMRestoreLocalVector(networkdm, &localXdot));
729 
730   PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F));
731   PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F));
732   PetscCall(DMRestoreLocalVector(networkdm, &localF));
733   PetscFunctionReturn(0);
734 }
735 
736 /* This function is used for solving the algebraic system only during fault on and
737    off times. It computes the entire F and then zeros out the part corresponding to
738    differential equations
739  F = [0;g(y)];
740 */
741 PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, void *ctx) {
742   DM                 networkdm;
743   Vec                localX, localF;
744   PetscInt           vfrom, vto, offsetfrom, offsetto;
745   PetscInt           v, vStart, vEnd, e;
746   PetscScalar       *farr;
747   Userctx           *user = (Userctx *)ctx;
748   const PetscScalar *xarr;
749   void              *component;
750   PetscScalar        Vr = 0, Vi = 0;
751 
752   PetscFunctionBegin;
753   PetscCall(VecSet(F, 0.0));
754   PetscCall(SNESGetDM(snes, &networkdm));
755   PetscCall(DMGetLocalVector(networkdm, &localF));
756   PetscCall(DMGetLocalVector(networkdm, &localX));
757   PetscCall(VecSet(localF, 0.0));
758 
759   /* update ghost values of locaX and locaXdot */
760   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
761   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
762 
763   PetscCall(VecGetArrayRead(localX, &xarr));
764   PetscCall(VecGetArray(localF, &farr));
765 
766   PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
767 
768   for (v = vStart; v < vEnd; v++) {
769     PetscInt    i, j, offsetbus, offsetgen, key, numComps;
770     PetscScalar Yffr, Yffi, Vm, Vm2, Vm0, Vr0 = 0, Vi0 = 0, PD, QD;
771     Bus        *bus;
772     Gen        *gen;
773     Load       *load;
774     PetscBool   ghostvtex;
775 
776     PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
777     PetscCall(DMNetworkGetNumComponents(networkdm, v, &numComps));
778 
779     for (j = 0; j < numComps; j++) {
780       PetscCall(DMNetworkGetComponent(networkdm, v, j, &key, &component, NULL));
781       if (key == 1) {
782         PetscInt        nconnedges;
783         const PetscInt *connedges;
784 
785         bus = (Bus *)(component);
786         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetbus));
787         if (!ghostvtex) {
788           Vr = xarr[offsetbus];
789           Vi = xarr[offsetbus + 1];
790 
791           Yffr = bus->yff[1];
792           Yffi = bus->yff[0];
793           if (user->alg_flg) {
794             Yffr += user->ybusfault[bus->id * 2 + 1];
795             Yffi += user->ybusfault[bus->id * 2];
796           }
797           Vr0 = bus->vr;
798           Vi0 = bus->vi;
799 
800           farr[offsetbus] += Yffi * Vr + Yffr * Vi;
801           farr[offsetbus + 1] += Yffr * Vr - Yffi * Vi;
802         }
803         PetscCall(DMNetworkGetSupportingEdges(networkdm, v, &nconnedges, &connedges));
804 
805         for (i = 0; i < nconnedges; i++) {
806           Branch         *branch;
807           PetscInt        keye;
808           PetscScalar     Yfti, Yftr, Vfr, Vfi, Vtr, Vti;
809           const PetscInt *cone;
810 
811           e = connedges[i];
812           PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, (void **)&branch, NULL));
813 
814           Yfti = branch->yft[0];
815           Yftr = branch->yft[1];
816 
817           PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone));
818           vfrom = cone[0];
819           vto   = cone[1];
820 
821           PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, 0, &offsetfrom));
822           PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, 0, &offsetto));
823 
824           /*From bus and to bus real and imaginary voltages */
825           Vfr = xarr[offsetfrom];
826           Vfi = xarr[offsetfrom + 1];
827           Vtr = xarr[offsetto];
828           Vti = xarr[offsetto + 1];
829 
830           if (vfrom == v) {
831             farr[offsetfrom] += Yftr * Vti + Yfti * Vtr;
832             farr[offsetfrom + 1] += Yftr * Vtr - Yfti * Vti;
833           } else {
834             farr[offsetto] += Yftr * Vfi + Yfti * Vfr;
835             farr[offsetto + 1] += Yftr * Vfr - Yfti * Vfi;
836           }
837         }
838       } else if (key == 2) {
839         if (!ghostvtex) {
840           PetscScalar Eqp, Edp, delta; /* Generator variables */
841           PetscScalar Id, Iq;          /* Generator dq axis currents */
842           PetscScalar Vd, Vq, IGr, IGi, Zdq_inv[4], det;
843           PetscScalar Xdp, Xqp, Rs; /* Generator parameters */
844 
845           gen = (Gen *)(component);
846           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetgen));
847 
848           /* Generator state variables */
849           Eqp   = xarr[offsetgen];
850           Edp   = xarr[offsetgen + 1];
851           delta = xarr[offsetgen + 2];
852           /* w     = xarr[idx+3]; not being used */
853           Id    = xarr[offsetgen + 4];
854           Iq    = xarr[offsetgen + 5];
855 
856           /* Generator parameters */
857           Xdp = gen->Xdp;
858           Xqp = gen->Xqp;
859           Rs  = gen->Rs;
860 
861           /* Set generator differential equation residual functions to zero */
862           farr[offsetgen]     = 0;
863           farr[offsetgen + 1] = 0;
864           farr[offsetgen + 2] = 0;
865           farr[offsetgen + 3] = 0;
866 
867           PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
868 
869           /* Algebraic equations for stator currents */
870           det = Rs * Rs + Xdp * Xqp;
871 
872           Zdq_inv[0] = Rs / det;
873           Zdq_inv[1] = Xqp / det;
874           Zdq_inv[2] = -Xdp / det;
875           Zdq_inv[3] = Rs / det;
876 
877           farr[offsetgen + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
878           farr[offsetgen + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
879 
880           /* Add generator current injection to network */
881           PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
882 
883           farr[offsetbus] -= IGi;
884           farr[offsetbus + 1] -= IGr;
885 
886           /* Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);*/ /* a compiler warning: "Value stored to 'Vm' is never read" - comment out by Hong Zhang */
887         }
888       } else if (key == 3) {
889         if (!ghostvtex) {
890           PetscInt offsetexc;
891           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetexc));
892           /* Set exciter differential equation residual functions equal to zero*/
893           farr[offsetexc]     = 0;
894           farr[offsetexc + 1] = 0;
895           farr[offsetexc + 2] = 0;
896         }
897       } else if (key == 4) {
898         if (!ghostvtex) {
899           PetscInt     k, ld_nsegsp, ld_nsegsq;
900           PetscScalar *ld_alphap, *ld_betap, *ld_alphaq, *ld_betaq, PD0, QD0, IDr, IDi;
901 
902           load = (Load *)(component);
903 
904           /* Load Parameters */
905           ld_nsegsp = load->ld_nsegsp;
906           ld_alphap = load->ld_alphap;
907           ld_betap  = load->ld_betap;
908           ld_nsegsq = load->ld_nsegsq;
909           ld_alphaq = load->ld_alphaq;
910           ld_betaq  = load->ld_betaq;
911 
912           PD0 = load->PD0;
913           QD0 = load->QD0;
914 
915           Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
916           Vm2 = Vm * Vm;
917           Vm0 = PetscSqrtScalar(Vr0 * Vr0 + Vi0 * Vi0);
918           PD = QD = 0.0;
919           for (k = 0; k < ld_nsegsp; k++) PD += ld_alphap[k] * PD0 * PetscPowScalar((Vm / Vm0), ld_betap[k]);
920           for (k = 0; k < ld_nsegsq; k++) QD += ld_alphaq[k] * QD0 * PetscPowScalar((Vm / Vm0), ld_betaq[k]);
921 
922           /* Load currents */
923           IDr = (PD * Vr + QD * Vi) / Vm2;
924           IDi = (-QD * Vr + PD * Vi) / Vm2;
925 
926           farr[offsetbus] += IDi;
927           farr[offsetbus + 1] += IDr;
928         }
929       }
930     }
931   }
932 
933   PetscCall(VecRestoreArrayRead(localX, &xarr));
934   PetscCall(VecRestoreArray(localF, &farr));
935   PetscCall(DMRestoreLocalVector(networkdm, &localX));
936 
937   PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F));
938   PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F));
939   PetscCall(DMRestoreLocalVector(networkdm, &localF));
940   PetscFunctionReturn(0);
941 }
942 
943 int main(int argc, char **argv) {
944   PetscInt    i, j, *edgelist = NULL, eStart, eEnd, vStart, vEnd;
945   PetscInt    genj, excj, loadj, componentkey[5];
946   PetscInt    nc = 1; /* No. of copies (default = 1) */
947   PetscMPIInt size, rank;
948   Vec         X, F_alg;
949   TS          ts;
950   SNES        snes_alg, snes;
951   Bus        *bus;
952   Branch     *branch;
953   Gen        *gen;
954   Exc        *exc;
955   Load       *load;
956   DM          networkdm;
957 #if defined(PETSC_USE_LOG)
958   PetscLogStage stage1;
959 #endif
960   Userctx  user;
961   KSP      ksp;
962   PC       pc;
963   PetscInt numEdges = 0;
964 
965   PetscFunctionBeginUser;
966   PetscCall(PetscInitialize(&argc, &argv, "ex9busnetworkops", help));
967   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nc", &nc, NULL));
968   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
969   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
970 
971   /* Read initial voltage vector and Ybus */
972   if (rank == 0) { PetscCall(read_data(nc, &gen, &exc, &load, &bus, &branch, &edgelist)); }
973 
974   PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &networkdm));
975   PetscCall(DMNetworkRegisterComponent(networkdm, "branchstruct", sizeof(Branch), &componentkey[0]));
976   PetscCall(DMNetworkRegisterComponent(networkdm, "busstruct", sizeof(Bus), &componentkey[1]));
977   PetscCall(DMNetworkRegisterComponent(networkdm, "genstruct", sizeof(Gen), &componentkey[2]));
978   PetscCall(DMNetworkRegisterComponent(networkdm, "excstruct", sizeof(Exc), &componentkey[3]));
979   PetscCall(DMNetworkRegisterComponent(networkdm, "loadstruct", sizeof(Load), &componentkey[4]));
980 
981   PetscCall(PetscLogStageRegister("Create network", &stage1));
982   PetscCall(PetscLogStagePush(stage1));
983 
984   /* Set local number of edges and edge connectivity */
985   if (rank == 0) numEdges = NBRANCH * nc + (nc - 1);
986   PetscCall(DMNetworkSetNumSubNetworks(networkdm, PETSC_DECIDE, 1));
987   PetscCall(DMNetworkAddSubnetwork(networkdm, NULL, numEdges, edgelist, NULL));
988 
989   /* Set up the network layout */
990   PetscCall(DMNetworkLayoutSetUp(networkdm));
991 
992   if (rank == 0) { PetscCall(PetscFree(edgelist)); }
993 
994   /* Add network components (physical parameters of nodes and branches) and number of variables */
995   if (rank == 0) {
996     PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd));
997     genj  = 0;
998     loadj = 0;
999     excj  = 0;
1000     for (i = eStart; i < eEnd; i++) { PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[0], &branch[i - eStart], 0)); }
1001 
1002     PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
1003 
1004     for (i = vStart; i < vEnd; i++) {
1005       PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[1], &bus[i - vStart], 2));
1006       if (bus[i - vStart].nofgen) {
1007         for (j = 0; j < bus[i - vStart].nofgen; j++) {
1008           /* Add generator */
1009           PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[2], &gen[genj++], 6));
1010           /* Add exciter */
1011           PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[3], &exc[excj++], 3));
1012         }
1013       }
1014       if (bus[i - vStart].nofload) {
1015         for (j = 0; j < bus[i - vStart].nofload; j++) { PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[4], &load[loadj++], 0)); }
1016       }
1017     }
1018   }
1019 
1020   PetscCall(DMSetUp(networkdm));
1021 
1022   if (rank == 0) { PetscCall(PetscFree5(bus, gen, load, branch, exc)); }
1023 
1024   /* for parallel options: Network partitioning and distribution of data */
1025   if (size > 1) { PetscCall(DMNetworkDistribute(&networkdm, 0)); }
1026   PetscCall(PetscLogStagePop());
1027 
1028   PetscCall(DMCreateGlobalVector(networkdm, &X));
1029 
1030   PetscCall(SetInitialGuess(networkdm, X));
1031 
1032   /* Options for fault simulation */
1033   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
1034   user.tfaulton  = 0.02;
1035   user.tfaultoff = 0.05;
1036   user.Rfault    = 0.0001;
1037   user.faultbus  = 8;
1038   PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
1039   PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
1040   PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
1041   user.t0   = 0.0;
1042   user.tmax = 0.1;
1043   PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
1044   PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
1045 
1046   PetscCall(PetscMalloc1(18 * nc, &user.ybusfault));
1047   for (i = 0; i < 18 * nc; i++) { user.ybusfault[i] = 0; }
1048   user.ybusfault[user.faultbus * 2 + 1] = 1 / user.Rfault;
1049   PetscOptionsEnd();
1050 
1051   /* Setup TS solver                                           */
1052   /*--------------------------------------------------------*/
1053   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
1054   PetscCall(TSSetDM(ts, (DM)networkdm));
1055   PetscCall(TSSetType(ts, TSCN));
1056 
1057   PetscCall(TSGetSNES(ts, &snes));
1058   PetscCall(SNESGetKSP(snes, &ksp));
1059   PetscCall(KSPGetPC(ksp, &pc));
1060   PetscCall(PCSetType(pc, PCBJACOBI));
1061 
1062   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)FormIFunction, &user));
1063   PetscCall(TSSetMaxTime(ts, user.tfaulton));
1064   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1065   PetscCall(TSSetTimeStep(ts, 0.01));
1066   PetscCall(TSSetFromOptions(ts));
1067 
1068   /*user.alg_flg = PETSC_TRUE is the period when fault exists. We add fault admittance to Ybus matrix.
1069     eg, fault bus is 8. Y88(new)=Y88(old)+Yfault. */
1070   user.alg_flg = PETSC_FALSE;
1071 
1072   /* Prefault period */
1073   if (rank == 0) { PetscCall(PetscPrintf(PETSC_COMM_SELF, "... (1) Prefault period ... \n")); }
1074 
1075   PetscCall(TSSetSolution(ts, X));
1076   PetscCall(TSSetUp(ts));
1077   PetscCall(TSSolve(ts, X));
1078 
1079   /* Create the nonlinear solver for solving the algebraic system */
1080   PetscCall(VecDuplicate(X, &F_alg));
1081 
1082   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
1083   PetscCall(SNESSetDM(snes_alg, (DM)networkdm));
1084   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, &user));
1085   PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_"));
1086   PetscCall(SNESSetFromOptions(snes_alg));
1087 
1088   /* Apply disturbance - resistive fault at user.faultbus */
1089   /* This is done by adding shunt conductance to the diagonal location
1090      in the Ybus matrix */
1091   user.alg_flg = PETSC_TRUE;
1092 
1093   /* Solve the algebraic equations */
1094   if (rank == 0) { PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n... (2) Apply disturbance, solve algebraic equations ... \n")); }
1095   PetscCall(SNESSolve(snes_alg, NULL, X));
1096 
1097   /* Disturbance period */
1098   PetscCall(TSSetTime(ts, user.tfaulton));
1099   PetscCall(TSSetMaxTime(ts, user.tfaultoff));
1100   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1101   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)FormIFunction, &user));
1102 
1103   user.alg_flg = PETSC_TRUE;
1104   if (rank == 0) { PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n... (3) Disturbance period ... \n")); }
1105   PetscCall(TSSolve(ts, X));
1106 
1107   /* Remove the fault */
1108   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, &user));
1109 
1110   user.alg_flg = PETSC_FALSE;
1111   /* Solve the algebraic equations */
1112   if (rank == 0) { PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n... (4) Remove fault, solve algebraic equations ... \n")); }
1113   PetscCall(SNESSolve(snes_alg, NULL, X));
1114   PetscCall(SNESDestroy(&snes_alg));
1115 
1116   /* Post-disturbance period */
1117   PetscCall(TSSetTime(ts, user.tfaultoff));
1118   PetscCall(TSSetMaxTime(ts, user.tmax));
1119   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
1120   PetscCall(TSSetIFunction(ts, NULL, (TSIFunction)FormIFunction, &user));
1121 
1122   user.alg_flg = PETSC_FALSE;
1123   if (rank == 0) { PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n... (5) Post-disturbance period ... \n")); }
1124   PetscCall(TSSolve(ts, X));
1125 
1126   PetscCall(PetscFree(user.ybusfault));
1127   PetscCall(VecDestroy(&F_alg));
1128   PetscCall(VecDestroy(&X));
1129   PetscCall(DMDestroy(&networkdm));
1130   PetscCall(TSDestroy(&ts));
1131   PetscCall(PetscFinalize());
1132   return 0;
1133 }
1134 
1135 /*TEST
1136 
1137    build:
1138       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
1139 
1140    test:
1141       args: -ts_monitor -snes_converged_reason -alg_snes_converged_reason
1142       localrunfiles: X.bin Ybus.bin ex9busnetworkops
1143 
1144    test:
1145       suffix: 2
1146       nsize: 2
1147       args: -ts_monitor -snes_converged_reason -alg_snes_converged_reason
1148       localrunfiles: X.bin Ybus.bin ex9busnetworkops
1149 
1150 TEST*/
1151