xref: /petsc/src/ts/tutorials/power_grid/stability_9bus/ex9busdmnetwork.c (revision 4e8208cbcbc709572b8abe32f33c78b69c819375)
1c4762a1bSJed Brown static char help[] = "This example uses the same problem set up of ex9busdmnetwork.c. \n\
2c4762a1bSJed Brown It demonstrates setting and accessing of variables for individual components, instead of \n\
3c4762a1bSJed Brown the network vertices (as used in ex9busdmnetwork.c). This is especially useful where vertices \n\
4c4762a1bSJed Brown /edges have multiple-components associated with them and one or more components has physics \n\
5c4762a1bSJed Brown associated with it. \n\
6c4762a1bSJed Brown Input parameters include:\n\
7c4762a1bSJed Brown   -nc : number of copies of the base case\n\n";
8c4762a1bSJed Brown 
91bb3edfdSPatrick Sanan /*
10c4762a1bSJed Brown    This example was modified from ex9busdmnetwork.c.
11c4762a1bSJed Brown */
12c4762a1bSJed Brown 
13c4762a1bSJed Brown #include <petscts.h>
14c4762a1bSJed Brown #include <petscdmnetwork.h>
15c4762a1bSJed Brown 
16c4762a1bSJed Brown #define FREQ    60
17c4762a1bSJed Brown #define W_S     (2 * PETSC_PI * FREQ)
18c4762a1bSJed Brown #define NGEN    3 /* No. of generators in the 9 bus system */
19c4762a1bSJed Brown #define NLOAD   3 /* No. of loads in the 9 bus system */
20c4762a1bSJed Brown #define NBUS    9 /* No. of buses in the 9 bus system */
21c4762a1bSJed Brown #define NBRANCH 9 /* No. of branches in the 9 bus system */
22c4762a1bSJed Brown 
23c4762a1bSJed Brown typedef struct {
24c4762a1bSJed Brown   PetscInt    id;    /* Bus Number or extended bus name*/
25c4762a1bSJed Brown   PetscScalar mbase; /* MVA base of the machine */
26c4762a1bSJed Brown   PetscScalar PG;    /* Generator active power output */
27c4762a1bSJed Brown   PetscScalar QG;    /* Generator reactive power output */
28c4762a1bSJed Brown 
29c4762a1bSJed Brown   /* Generator constants */
30c4762a1bSJed Brown   PetscScalar H;    /* Inertia constant */
31c4762a1bSJed Brown   PetscScalar Rs;   /* Stator Resistance */
32c4762a1bSJed Brown   PetscScalar Xd;   /* d-axis reactance */
33c4762a1bSJed Brown   PetscScalar Xdp;  /* d-axis transient reactance */
34c4762a1bSJed Brown   PetscScalar Xq;   /* q-axis reactance Xq(1) set to 0.4360, value given in text 0.0969 */
35c4762a1bSJed Brown   PetscScalar Xqp;  /* q-axis transient reactance */
36c4762a1bSJed Brown   PetscScalar Td0p; /* d-axis open circuit time constant */
37c4762a1bSJed Brown   PetscScalar Tq0p; /* q-axis open circuit time constant */
38c4762a1bSJed Brown   PetscScalar M;    /* M = 2*H/W_S */
39c4762a1bSJed Brown   PetscScalar D;    /* D = 0.1*M */
40c4762a1bSJed Brown   PetscScalar TM;   /* Mechanical Torque */
41c4762a1bSJed Brown } Gen;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown typedef struct {
44c4762a1bSJed Brown   /* Exciter system constants */
45c4762a1bSJed Brown   PetscScalar KA;     /* Voltage regulator gain constant */
46c4762a1bSJed Brown   PetscScalar TA;     /* Voltage regulator time constant */
47c4762a1bSJed Brown   PetscScalar KE;     /* Exciter gain constant */
48c4762a1bSJed Brown   PetscScalar TE;     /* Exciter time constant */
49c4762a1bSJed Brown   PetscScalar KF;     /* Feedback stabilizer gain constant */
50c4762a1bSJed Brown   PetscScalar TF;     /* Feedback stabilizer time constant */
51c4762a1bSJed Brown   PetscScalar k1, k2; /* calculating the saturation function SE = k1*exp(k2*Efd) */
52c4762a1bSJed Brown   PetscScalar Vref;   /* Voltage regulator voltage setpoint */
53c4762a1bSJed Brown } Exc;
54c4762a1bSJed Brown 
55c4762a1bSJed Brown typedef struct {
56c4762a1bSJed Brown   PetscInt    id;      /* node id */
57c4762a1bSJed Brown   PetscInt    nofgen;  /* Number of generators at the bus*/
58c4762a1bSJed Brown   PetscInt    nofload; /*  Number of load at the bus*/
59c4762a1bSJed Brown   PetscScalar yff[2];  /* yff[0]= imaginary part of admittance, yff[1]=real part of admittance*/
60c4762a1bSJed Brown   PetscScalar vr;      /* Real component of bus voltage */
61c4762a1bSJed Brown   PetscScalar vi;      /* Imaginary component of bus voltage */
62c4762a1bSJed Brown } Bus;
63c4762a1bSJed Brown 
64c4762a1bSJed Brown /* Load constants
65c4762a1bSJed Brown   We use a composite load model that describes the load and reactive powers at each time instant as follows
66c4762a1bSJed Brown   P(t) = \sum\limits_{i=0}^ld_nsegsp \ld_alphap_i*P_D0(\frac{V_m(t)}{V_m0})^\ld_betap_i
67c4762a1bSJed Brown   Q(t) = \sum\limits_{i=0}^ld_nsegsq \ld_alphaq_i*Q_D0(\frac{V_m(t)}{V_m0})^\ld_betaq_i
68c4762a1bSJed Brown   where
69c4762a1bSJed Brown     id                  - index of the load
70c4762a1bSJed Brown     ld_nsegsp,ld_nsegsq - Number of individual load models for real and reactive power loads
71c4762a1bSJed Brown     ld_alphap,ld_alphap - Percentage contribution (weights) or loads
72c4762a1bSJed Brown     P_D0                - Real power load
73c4762a1bSJed Brown     Q_D0                - Reactive power load
74c4762a1bSJed Brown     Vm(t)              - Voltage magnitude at time t
75c4762a1bSJed Brown     Vm0                - Voltage magnitude at t = 0
76c4762a1bSJed Brown     ld_betap, ld_betaq  - exponents describing the load model for real and reactive part
77c4762a1bSJed Brown 
78c4762a1bSJed Brown     Note: All loads have the same characteristic currently.
79c4762a1bSJed Brown   */
80c4762a1bSJed Brown typedef struct {
81c4762a1bSJed Brown   PetscInt    id; /* bus id */
82c4762a1bSJed Brown   PetscInt    ld_nsegsp, ld_nsegsq;
83c4762a1bSJed Brown   PetscScalar PD0, QD0;
84c4762a1bSJed Brown   PetscScalar ld_alphap[3]; /* ld_alphap=[1,0,0], an array, not a value, so use *ld_alphap; */
85c4762a1bSJed Brown   PetscScalar ld_betap[3], ld_alphaq[3], ld_betaq[3];
86c4762a1bSJed Brown } Load;
87c4762a1bSJed Brown 
88c4762a1bSJed Brown typedef struct {
89c4762a1bSJed Brown   PetscInt    id;     /* node id */
90c4762a1bSJed Brown   PetscScalar yft[2]; /* yft[0]= imaginary part of admittance, yft[1]=real part of admittance*/
91c4762a1bSJed Brown } Branch;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown typedef struct {
94c4762a1bSJed Brown   PetscReal    tfaulton, tfaultoff; /* Fault on and off times */
95c4762a1bSJed Brown   PetscReal    t;
96c4762a1bSJed Brown   PetscReal    t0, tmax; /* initial time and final time */
97c4762a1bSJed Brown   PetscInt     faultbus; /* Fault bus */
98c4762a1bSJed Brown   PetscScalar  Rfault;   /* Fault resistance (pu) */
99c4762a1bSJed Brown   PetscScalar *ybusfault;
100c4762a1bSJed Brown   PetscBool    alg_flg;
101c4762a1bSJed Brown } Userctx;
102c4762a1bSJed Brown 
103c4762a1bSJed Brown /* Used to read data into the DMNetwork components */
read_data(PetscInt nc,Gen ** pgen,Exc ** pexc,Load ** pload,Bus ** pbus,Branch ** pbranch,PetscInt ** pedgelist)104d71ae5a4SJacob Faibussowitsch PetscErrorCode read_data(PetscInt nc, Gen **pgen, Exc **pexc, Load **pload, Bus **pbus, Branch **pbranch, PetscInt **pedgelist)
105d71ae5a4SJacob Faibussowitsch {
106c4762a1bSJed Brown   PetscInt           i, j, row[1], col[2];
107c4762a1bSJed Brown   PetscInt          *edgelist;
108c4762a1bSJed Brown   PetscInt           nofgen[9]  = {1, 1, 1, 0, 0, 0, 0, 0, 0}; /* Buses at which generators are incident */
109c4762a1bSJed Brown   PetscInt           nofload[9] = {0, 0, 0, 0, 1, 1, 0, 1, 0}; /* Buses at which loads are incident */
110c4762a1bSJed Brown   const PetscScalar *varr;
111c4762a1bSJed Brown   PetscScalar        M[3], D[3];
112c4762a1bSJed Brown   Bus               *bus;
113c4762a1bSJed Brown   Branch            *branch;
114c4762a1bSJed Brown   Gen               *gen;
115c4762a1bSJed Brown   Exc               *exc;
116c4762a1bSJed Brown   Load              *load;
117c4762a1bSJed Brown   Mat                Ybus;
118c4762a1bSJed Brown   Vec                V0;
119c4762a1bSJed Brown 
120c4762a1bSJed Brown   /*10 parameters*/
121c4762a1bSJed Brown   /* Generator real and reactive powers (found via loadflow) */
122c4762a1bSJed Brown   static const PetscScalar PG[3] = {0.716786142395021, 1.630000000000000, 0.850000000000000};
123c4762a1bSJed Brown   static const PetscScalar QG[3] = {0.270702180178785, 0.066120127797275, -0.108402221791588};
124c4762a1bSJed Brown 
125c4762a1bSJed Brown   /* Generator constants */
126c4762a1bSJed Brown   static const PetscScalar H[3]    = {23.64, 6.4, 3.01};       /* Inertia constant */
127c4762a1bSJed Brown   static const PetscScalar Rs[3]   = {0.0, 0.0, 0.0};          /* Stator Resistance */
128c4762a1bSJed Brown   static const PetscScalar Xd[3]   = {0.146, 0.8958, 1.3125};  /* d-axis reactance */
129c4762a1bSJed Brown   static const PetscScalar Xdp[3]  = {0.0608, 0.1198, 0.1813}; /* d-axis transient reactance */
130c4762a1bSJed Brown   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 */
131c4762a1bSJed Brown   static const PetscScalar Xqp[3]  = {0.0969, 0.1969, 0.25};   /* q-axis transient reactance */
132c4762a1bSJed Brown   static const PetscScalar Td0p[3] = {8.96, 6.0, 5.89};        /* d-axis open circuit time constant */
133c4762a1bSJed Brown   static const PetscScalar Tq0p[3] = {0.31, 0.535, 0.6};       /* q-axis open circuit time constant */
134c4762a1bSJed Brown 
135c4762a1bSJed Brown   /* Exciter system constants (8 parameters)*/
136c4762a1bSJed Brown   static const PetscScalar KA[3] = {20.0, 20.0, 20.0};    /* Voltage regulartor gain constant */
137c4762a1bSJed Brown   static const PetscScalar TA[3] = {0.2, 0.2, 0.2};       /* Voltage regulator time constant */
138c4762a1bSJed Brown   static const PetscScalar KE[3] = {1.0, 1.0, 1.0};       /* Exciter gain constant */
139c4762a1bSJed Brown   static const PetscScalar TE[3] = {0.314, 0.314, 0.314}; /* Exciter time constant */
140c4762a1bSJed Brown   static const PetscScalar KF[3] = {0.063, 0.063, 0.063}; /* Feedback stabilizer gain constant */
141c4762a1bSJed Brown   static const PetscScalar TF[3] = {0.35, 0.35, 0.35};    /* Feedback stabilizer time constant */
142c4762a1bSJed Brown   static const PetscScalar k1[3] = {0.0039, 0.0039, 0.0039};
143c4762a1bSJed Brown   static const PetscScalar k2[3] = {1.555, 1.555, 1.555}; /* k1 and k2 for calculating the saturation function SE = k1*exp(k2*Efd) */
144c4762a1bSJed Brown 
145c4762a1bSJed Brown   /* Load constants */
146c4762a1bSJed Brown   static const PetscScalar PD0[3]       = {1.25, 0.9, 1.0};
147c4762a1bSJed Brown   static const PetscScalar QD0[3]       = {0.5, 0.3, 0.35};
148c4762a1bSJed Brown   static const PetscScalar ld_alphaq[3] = {1, 0, 0};
149c4762a1bSJed Brown   static const PetscScalar ld_betaq[3]  = {2, 1, 0};
150c4762a1bSJed Brown   static const PetscScalar ld_betap[3]  = {2, 1, 0};
151c4762a1bSJed Brown   static const PetscScalar ld_alphap[3] = {1, 0, 0};
152c4762a1bSJed Brown   PetscInt                 ld_nsegsp[3] = {3, 3, 3};
153c4762a1bSJed Brown   PetscInt                 ld_nsegsq[3] = {3, 3, 3};
154c4762a1bSJed Brown   PetscViewer              Xview, Ybusview;
155c4762a1bSJed Brown   PetscInt                 neqs_net, m, n;
156c4762a1bSJed Brown 
157c4762a1bSJed Brown   PetscFunctionBeginUser;
158c4762a1bSJed Brown   /* Read V0 and Ybus from files */
1599566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "X.bin", FILE_MODE_READ, &Xview));
1609566063dSJacob Faibussowitsch   PetscCall(PetscViewerBinaryOpen(PETSC_COMM_SELF, "Ybus.bin", FILE_MODE_READ, &Ybusview));
1619566063dSJacob Faibussowitsch   PetscCall(VecCreate(PETSC_COMM_SELF, &V0));
1629566063dSJacob Faibussowitsch   PetscCall(VecLoad(V0, Xview));
163c4762a1bSJed Brown 
1649566063dSJacob Faibussowitsch   PetscCall(MatCreate(PETSC_COMM_SELF, &Ybus));
1659566063dSJacob Faibussowitsch   PetscCall(MatSetType(Ybus, MATBAIJ));
1669566063dSJacob Faibussowitsch   PetscCall(MatLoad(Ybus, Ybusview));
167c4762a1bSJed Brown 
168c4762a1bSJed Brown   /* Destroy unnecessary stuff */
1699566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&Xview));
1709566063dSJacob Faibussowitsch   PetscCall(PetscViewerDestroy(&Ybusview));
171c4762a1bSJed Brown 
1729566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(Ybus, &m, &n));
173c4762a1bSJed Brown   neqs_net = 2 * NBUS; /* # eqs. for network subsystem   */
174cad9d221SBarry Smith   PetscCheck(m == neqs_net && n == neqs_net, PETSC_COMM_SELF, PETSC_ERR_ARG_SIZ, "matrix Ybus is in wrong sizes");
175c4762a1bSJed Brown 
176c4762a1bSJed Brown   M[0] = 2 * H[0] / W_S;
177c4762a1bSJed Brown   M[1] = 2 * H[1] / W_S;
178c4762a1bSJed Brown   M[2] = 2 * H[2] / W_S;
179c4762a1bSJed Brown   D[0] = 0.1 * M[0];
180c4762a1bSJed Brown   D[1] = 0.1 * M[1];
181c4762a1bSJed Brown   D[2] = 0.1 * M[2];
182c4762a1bSJed Brown 
183d5b43468SJose E. Roman   /* Allocate memory for bus, generators, exciter, loads and branches */
1849566063dSJacob Faibussowitsch   PetscCall(PetscCalloc5(NBUS * nc, &bus, NGEN * nc, &gen, NLOAD * nc, &load, NBRANCH * nc + (nc - 1), &branch, NGEN * nc, &exc));
185c4762a1bSJed Brown 
1869566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(V0, &varr));
187c4762a1bSJed Brown 
188c4762a1bSJed Brown   /* read bus data */
189c4762a1bSJed Brown   for (i = 0; i < nc; i++) {
190c4762a1bSJed Brown     for (j = 0; j < NBUS; j++) {
191c4762a1bSJed Brown       bus[i * 9 + j].id      = i * 9 + j;
192c4762a1bSJed Brown       bus[i * 9 + j].nofgen  = nofgen[j];
193c4762a1bSJed Brown       bus[i * 9 + j].nofload = nofload[j];
194c4762a1bSJed Brown       bus[i * 9 + j].vr      = varr[2 * j];
195c4762a1bSJed Brown       bus[i * 9 + j].vi      = varr[2 * j + 1];
196c4762a1bSJed Brown       row[0]                 = 2 * j;
197c4762a1bSJed Brown       col[0]                 = 2 * j;
198c4762a1bSJed Brown       col[1]                 = 2 * j + 1;
199c4762a1bSJed Brown       /* real and imaginary part of admittance from Ybus into yff */
2009566063dSJacob Faibussowitsch       PetscCall(MatGetValues(Ybus, 1, row, 2, col, bus[i * 9 + j].yff));
201c4762a1bSJed Brown     }
202c4762a1bSJed Brown   }
203c4762a1bSJed Brown 
204c4762a1bSJed Brown   /* read generator data */
205c4762a1bSJed Brown   for (i = 0; i < nc; i++) {
206c4762a1bSJed Brown     for (j = 0; j < NGEN; j++) {
207c4762a1bSJed Brown       gen[i * 3 + j].id   = i * 3 + j;
208c4762a1bSJed Brown       gen[i * 3 + j].PG   = PG[j];
209c4762a1bSJed Brown       gen[i * 3 + j].QG   = QG[j];
210c4762a1bSJed Brown       gen[i * 3 + j].H    = H[j];
211c4762a1bSJed Brown       gen[i * 3 + j].Rs   = Rs[j];
212c4762a1bSJed Brown       gen[i * 3 + j].Xd   = Xd[j];
213c4762a1bSJed Brown       gen[i * 3 + j].Xdp  = Xdp[j];
214c4762a1bSJed Brown       gen[i * 3 + j].Xq   = Xq[j];
215c4762a1bSJed Brown       gen[i * 3 + j].Xqp  = Xqp[j];
216c4762a1bSJed Brown       gen[i * 3 + j].Td0p = Td0p[j];
217c4762a1bSJed Brown       gen[i * 3 + j].Tq0p = Tq0p[j];
218c4762a1bSJed Brown       gen[i * 3 + j].M    = M[j];
219c4762a1bSJed Brown       gen[i * 3 + j].D    = D[j];
220c4762a1bSJed Brown     }
221c4762a1bSJed Brown   }
222c4762a1bSJed Brown 
223c4762a1bSJed Brown   for (i = 0; i < nc; i++) {
224c4762a1bSJed Brown     for (j = 0; j < NGEN; j++) {
225c4762a1bSJed Brown       /* exciter system */
226c4762a1bSJed Brown       exc[i * 3 + j].KA = KA[j];
227c4762a1bSJed Brown       exc[i * 3 + j].TA = TA[j];
228c4762a1bSJed Brown       exc[i * 3 + j].KE = KE[j];
229c4762a1bSJed Brown       exc[i * 3 + j].TE = TE[j];
230c4762a1bSJed Brown       exc[i * 3 + j].KF = KF[j];
231c4762a1bSJed Brown       exc[i * 3 + j].TF = TF[j];
232c4762a1bSJed Brown       exc[i * 3 + j].k1 = k1[j];
233c4762a1bSJed Brown       exc[i * 3 + j].k2 = k2[j];
234c4762a1bSJed Brown     }
235c4762a1bSJed Brown   }
236c4762a1bSJed Brown 
237c4762a1bSJed Brown   /* read load data */
238c4762a1bSJed Brown   for (i = 0; i < nc; i++) {
239c4762a1bSJed Brown     for (j = 0; j < NLOAD; j++) {
240c4762a1bSJed Brown       load[i * 3 + j].id        = i * 3 + j;
241c4762a1bSJed Brown       load[i * 3 + j].PD0       = PD0[j];
242c4762a1bSJed Brown       load[i * 3 + j].QD0       = QD0[j];
243c4762a1bSJed Brown       load[i * 3 + j].ld_nsegsp = ld_nsegsp[j];
244c4762a1bSJed Brown 
245c4762a1bSJed Brown       load[i * 3 + j].ld_alphap[0] = ld_alphap[0];
246c4762a1bSJed Brown       load[i * 3 + j].ld_alphap[1] = ld_alphap[1];
247c4762a1bSJed Brown       load[i * 3 + j].ld_alphap[2] = ld_alphap[2];
248c4762a1bSJed Brown 
249c4762a1bSJed Brown       load[i * 3 + j].ld_alphaq[0] = ld_alphaq[0];
250c4762a1bSJed Brown       load[i * 3 + j].ld_alphaq[1] = ld_alphaq[1];
251c4762a1bSJed Brown       load[i * 3 + j].ld_alphaq[2] = ld_alphaq[2];
252c4762a1bSJed Brown 
253c4762a1bSJed Brown       load[i * 3 + j].ld_betap[0] = ld_betap[0];
254c4762a1bSJed Brown       load[i * 3 + j].ld_betap[1] = ld_betap[1];
255c4762a1bSJed Brown       load[i * 3 + j].ld_betap[2] = ld_betap[2];
256c4762a1bSJed Brown       load[i * 3 + j].ld_nsegsq   = ld_nsegsq[j];
257c4762a1bSJed Brown 
258c4762a1bSJed Brown       load[i * 3 + j].ld_betaq[0] = ld_betaq[0];
259c4762a1bSJed Brown       load[i * 3 + j].ld_betaq[1] = ld_betaq[1];
260c4762a1bSJed Brown       load[i * 3 + j].ld_betaq[2] = ld_betaq[2];
261c4762a1bSJed Brown     }
262c4762a1bSJed Brown   }
2639566063dSJacob Faibussowitsch   PetscCall(PetscCalloc1(2 * NBRANCH * nc + 2 * (nc - 1), &edgelist));
264c4762a1bSJed Brown 
265c4762a1bSJed Brown   /* read edgelist */
266c4762a1bSJed Brown   for (i = 0; i < nc; i++) {
267c4762a1bSJed Brown     for (j = 0; j < NBRANCH; j++) {
268c4762a1bSJed Brown       switch (j) {
269c4762a1bSJed Brown       case 0:
270c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 0 + 9 * i;
271c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 3 + 9 * i;
272c4762a1bSJed Brown         break;
273c4762a1bSJed Brown       case 1:
274c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 1 + 9 * i;
275c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 6 + 9 * i;
276c4762a1bSJed Brown         break;
277c4762a1bSJed Brown       case 2:
278c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 2 + 9 * i;
279c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 8 + 9 * i;
280c4762a1bSJed Brown         break;
281c4762a1bSJed Brown       case 3:
282c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 3 + 9 * i;
283c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 4 + 9 * i;
284c4762a1bSJed Brown         break;
285c4762a1bSJed Brown       case 4:
286c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 3 + 9 * i;
287c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 5 + 9 * i;
288c4762a1bSJed Brown         break;
289c4762a1bSJed Brown       case 5:
290c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 4 + 9 * i;
291c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 6 + 9 * i;
292c4762a1bSJed Brown         break;
293c4762a1bSJed Brown       case 6:
294c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 5 + 9 * i;
295c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 8 + 9 * i;
296c4762a1bSJed Brown         break;
297c4762a1bSJed Brown       case 7:
298c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 6 + 9 * i;
299c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 7 + 9 * i;
300c4762a1bSJed Brown         break;
301c4762a1bSJed Brown       case 8:
302c4762a1bSJed Brown         edgelist[i * 18 + 2 * j]     = 7 + 9 * i;
303c4762a1bSJed Brown         edgelist[i * 18 + 2 * j + 1] = 8 + 9 * i;
304c4762a1bSJed Brown         break;
305d71ae5a4SJacob Faibussowitsch       default:
306d71ae5a4SJacob Faibussowitsch         break;
307c4762a1bSJed Brown       }
308c4762a1bSJed Brown     }
309c4762a1bSJed Brown   }
310c4762a1bSJed Brown 
311c4762a1bSJed Brown   /* 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 */
312c4762a1bSJed Brown   for (i = 1; i < nc; i++) {
313c4762a1bSJed Brown     edgelist[18 * nc + 2 * (i - 1)]     = 8 + (i - 1) * 9;
314c4762a1bSJed Brown     edgelist[18 * nc + 2 * (i - 1) + 1] = 9 * i;
315c4762a1bSJed Brown 
316c4762a1bSJed Brown     /* adding admittances to the off-diagonal elements */
317c4762a1bSJed Brown     branch[9 * nc + (i - 1)].id     = 9 * nc + (i - 1);
318c4762a1bSJed Brown     branch[9 * nc + (i - 1)].yft[0] = 17.3611;
319c4762a1bSJed Brown     branch[9 * nc + (i - 1)].yft[1] = -0.0301407;
320c4762a1bSJed Brown 
321c4762a1bSJed Brown     /* subtracting admittances from the diagonal elements */
322c4762a1bSJed Brown     bus[9 * i - 1].yff[0] -= 17.3611;
323c4762a1bSJed Brown     bus[9 * i - 1].yff[1] -= -0.0301407;
324c4762a1bSJed Brown 
325c4762a1bSJed Brown     bus[9 * i].yff[0] -= 17.3611;
326c4762a1bSJed Brown     bus[9 * i].yff[1] -= -0.0301407;
327c4762a1bSJed Brown   }
328c4762a1bSJed Brown 
329c4762a1bSJed Brown   /* read branch data */
330c4762a1bSJed Brown   for (i = 0; i < nc; i++) {
331c4762a1bSJed Brown     for (j = 0; j < NBRANCH; j++) {
332c4762a1bSJed Brown       branch[i * 9 + j].id = i * 9 + j;
333c4762a1bSJed Brown 
334c4762a1bSJed Brown       row[0] = edgelist[2 * j] * 2;
335c4762a1bSJed Brown       col[0] = edgelist[2 * j + 1] * 2;
336c4762a1bSJed Brown       col[1] = edgelist[2 * j + 1] * 2 + 1;
3379566063dSJacob Faibussowitsch       PetscCall(MatGetValues(Ybus, 1, row, 2, col, branch[i * 9 + j].yft)); /*imaginary part of admittance*/
338c4762a1bSJed Brown     }
339c4762a1bSJed Brown   }
340c4762a1bSJed Brown 
341c4762a1bSJed Brown   *pgen      = gen;
342c4762a1bSJed Brown   *pexc      = exc;
343c4762a1bSJed Brown   *pload     = load;
344c4762a1bSJed Brown   *pbus      = bus;
345c4762a1bSJed Brown   *pbranch   = branch;
346c4762a1bSJed Brown   *pedgelist = edgelist;
347c4762a1bSJed Brown 
3489566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(V0, &varr));
349c4762a1bSJed Brown 
350c4762a1bSJed Brown   /* Destroy unnecessary stuff */
3519566063dSJacob Faibussowitsch   PetscCall(MatDestroy(&Ybus));
3529566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&V0));
3533ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
354c4762a1bSJed Brown }
355c4762a1bSJed Brown 
SetInitialGuess(DM networkdm,Vec X)356d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialGuess(DM networkdm, Vec X)
357d71ae5a4SJacob Faibussowitsch {
358c4762a1bSJed Brown   Bus         *bus;
359c4762a1bSJed Brown   Gen         *gen;
360c4762a1bSJed Brown   Exc         *exc;
361c4762a1bSJed Brown   PetscInt     v, vStart, vEnd, offset;
362c4762a1bSJed Brown   PetscInt     key, numComps, j;
363c4762a1bSJed Brown   PetscBool    ghostvtex;
364c4762a1bSJed Brown   Vec          localX;
365c4762a1bSJed Brown   PetscScalar *xarr;
366c4762a1bSJed Brown   PetscScalar  Vr = 0, Vi = 0, Vm = 0, Vm2; /* Terminal voltage variables */
367c4762a1bSJed Brown   PetscScalar  IGr, IGi;                    /* Generator real and imaginary current */
368c4762a1bSJed Brown   PetscScalar  Eqp, Edp, delta;             /* Generator variables */
369c4762a1bSJed Brown   PetscScalar  Efd = 0, RF, VR;             /* Exciter variables */
370c4762a1bSJed Brown   PetscScalar  Vd, Vq;                      /* Generator dq axis voltages */
371c4762a1bSJed Brown   PetscScalar  Id, Iq;                      /* Generator dq axis currents */
372c4762a1bSJed Brown   PetscScalar  theta;                       /* Generator phase angle */
373c4762a1bSJed Brown   PetscScalar  SE;
374c4762a1bSJed Brown   void        *component;
375c4762a1bSJed Brown 
376c4762a1bSJed Brown   PetscFunctionBegin;
3779566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
3789566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localX));
379c4762a1bSJed Brown 
3809566063dSJacob Faibussowitsch   PetscCall(VecSet(X, 0.0));
3819566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
3829566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
383c4762a1bSJed Brown 
3849566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX, &xarr));
385c4762a1bSJed Brown 
386c4762a1bSJed Brown   for (v = vStart; v < vEnd; v++) {
3879566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
388c4762a1bSJed Brown     if (ghostvtex) continue;
389c4762a1bSJed Brown 
3909566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm, v, &numComps));
391c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
3929566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(networkdm, v, j, &key, &component, NULL));
393c4762a1bSJed Brown       if (key == 1) {
394c4762a1bSJed Brown         bus = (Bus *)(component);
395c4762a1bSJed Brown 
3969566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offset));
397c4762a1bSJed Brown         xarr[offset]     = bus->vr;
398c4762a1bSJed Brown         xarr[offset + 1] = bus->vi;
399c4762a1bSJed Brown 
400c4762a1bSJed Brown         Vr = bus->vr;
401c4762a1bSJed Brown         Vi = bus->vi;
402c4762a1bSJed Brown       } else if (key == 2) {
403c4762a1bSJed Brown         gen = (Gen *)(component);
4049566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offset));
405c4762a1bSJed Brown         Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
406c4762a1bSJed Brown         Vm2 = Vm * Vm;
407c4762a1bSJed Brown         /* Real part of gen current */
408c4762a1bSJed Brown         IGr = (Vr * gen->PG + Vi * gen->QG) / Vm2;
409c4762a1bSJed Brown         /* Imaginary part of gen current */
410c4762a1bSJed Brown         IGi = (Vi * gen->PG - Vr * gen->QG) / Vm2;
411c4762a1bSJed Brown 
412c4762a1bSJed Brown         /* Machine angle */
413c4762a1bSJed Brown         delta = atan2(Vi + gen->Xq * IGr, Vr - gen->Xq * IGi);
414c4762a1bSJed Brown         theta = PETSC_PI / 2.0 - delta;
415c4762a1bSJed Brown 
416c4762a1bSJed Brown         /* d-axis stator current */
417c4762a1bSJed Brown         Id = IGr * PetscCosScalar(theta) - IGi * PetscSinScalar(theta);
418c4762a1bSJed Brown 
419c4762a1bSJed Brown         /* q-axis stator current */
420c4762a1bSJed Brown         Iq = IGr * PetscSinScalar(theta) + IGi * PetscCosScalar(theta);
421c4762a1bSJed Brown 
422c4762a1bSJed Brown         Vd = Vr * PetscCosScalar(theta) - Vi * PetscSinScalar(theta);
423c4762a1bSJed Brown         Vq = Vr * PetscSinScalar(theta) + Vi * PetscCosScalar(theta);
424c4762a1bSJed Brown 
425c4762a1bSJed Brown         /* d-axis transient EMF */
426c4762a1bSJed Brown         Edp = Vd + gen->Rs * Id - gen->Xqp * Iq;
427c4762a1bSJed Brown 
428c4762a1bSJed Brown         /* q-axis transient EMF */
429c4762a1bSJed Brown         Eqp = Vq + gen->Rs * Iq + gen->Xdp * Id;
430c4762a1bSJed Brown 
431c4762a1bSJed Brown         gen->TM = gen->PG;
432c4762a1bSJed Brown 
433c4762a1bSJed Brown         xarr[offset]     = Eqp;
434c4762a1bSJed Brown         xarr[offset + 1] = Edp;
435c4762a1bSJed Brown         xarr[offset + 2] = delta;
436c4762a1bSJed Brown         xarr[offset + 3] = W_S;
437c4762a1bSJed Brown         xarr[offset + 4] = Id;
438c4762a1bSJed Brown         xarr[offset + 5] = Iq;
439c4762a1bSJed Brown 
440c4762a1bSJed Brown         Efd = Eqp + (gen->Xd - gen->Xdp) * Id;
441c4762a1bSJed Brown 
442c4762a1bSJed Brown       } else if (key == 3) {
443c4762a1bSJed Brown         exc = (Exc *)(component);
4449566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offset));
445c4762a1bSJed Brown 
446c4762a1bSJed Brown         SE = exc->k1 * PetscExpScalar(exc->k2 * Efd);
447c4762a1bSJed Brown         VR = exc->KE * Efd + SE;
448c4762a1bSJed Brown         RF = exc->KF * Efd / exc->TF;
449c4762a1bSJed Brown 
450c4762a1bSJed Brown         xarr[offset]     = Efd;
451c4762a1bSJed Brown         xarr[offset + 1] = RF;
452c4762a1bSJed Brown         xarr[offset + 2] = VR;
453c4762a1bSJed Brown 
454c4762a1bSJed Brown         exc->Vref = Vm + (VR / exc->KA);
455c4762a1bSJed Brown       }
456c4762a1bSJed Brown     }
457c4762a1bSJed Brown   }
4589566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX, &xarr));
4599566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(networkdm, localX, ADD_VALUES, X));
4609566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(networkdm, localX, ADD_VALUES, X));
4619566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localX));
4623ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
463c4762a1bSJed Brown }
464c4762a1bSJed Brown 
465c4762a1bSJed Brown /* Converts from machine frame (dq) to network (phase a real,imag) reference frame */
dq2ri(PetscScalar Fd,PetscScalar Fq,PetscScalar delta,PetscScalar * Fr,PetscScalar * Fi)466d71ae5a4SJacob Faibussowitsch PetscErrorCode dq2ri(PetscScalar Fd, PetscScalar Fq, PetscScalar delta, PetscScalar *Fr, PetscScalar *Fi)
467d71ae5a4SJacob Faibussowitsch {
468c4762a1bSJed Brown   PetscFunctionBegin;
469c4762a1bSJed Brown   *Fr = Fd * PetscSinScalar(delta) + Fq * PetscCosScalar(delta);
470c4762a1bSJed Brown   *Fi = -Fd * PetscCosScalar(delta) + Fq * PetscSinScalar(delta);
4713ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
472c4762a1bSJed Brown }
473c4762a1bSJed Brown 
474c4762a1bSJed Brown /* Converts from network frame ([phase a real,imag) to machine (dq) reference frame */
ri2dq(PetscScalar Fr,PetscScalar Fi,PetscScalar delta,PetscScalar * Fd,PetscScalar * Fq)475d71ae5a4SJacob Faibussowitsch PetscErrorCode ri2dq(PetscScalar Fr, PetscScalar Fi, PetscScalar delta, PetscScalar *Fd, PetscScalar *Fq)
476d71ae5a4SJacob Faibussowitsch {
477c4762a1bSJed Brown   PetscFunctionBegin;
478c4762a1bSJed Brown   *Fd = Fr * PetscSinScalar(delta) - Fi * PetscCosScalar(delta);
479c4762a1bSJed Brown   *Fq = Fr * PetscCosScalar(delta) + Fi * PetscSinScalar(delta);
4803ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
481c4762a1bSJed Brown }
482c4762a1bSJed Brown 
483c4762a1bSJed Brown /* Computes F(t,U,U_t) where F() = 0 is the DAE to be solved. */
FormIFunction(TS ts,PetscReal t,Vec X,Vec Xdot,Vec F,Userctx * user)484d71ae5a4SJacob Faibussowitsch PetscErrorCode FormIFunction(TS ts, PetscReal t, Vec X, Vec Xdot, Vec F, Userctx *user)
485d71ae5a4SJacob Faibussowitsch {
486c4762a1bSJed Brown   DM                 networkdm;
487c4762a1bSJed Brown   Vec                localX, localXdot, localF;
488c4762a1bSJed Brown   PetscInt           vfrom, vto, offsetfrom, offsetto;
489c4762a1bSJed Brown   PetscInt           v, vStart, vEnd, e;
490c4762a1bSJed Brown   PetscScalar       *farr;
4917b5e0a9aSStefano Zampini   PetscScalar        Vd = 0, Vq = 0, SE;
492c4762a1bSJed Brown   const PetscScalar *xarr, *xdotarr;
493c4762a1bSJed Brown   void              *component;
494c4762a1bSJed Brown   PetscScalar        Vr = 0, Vi = 0;
495c4762a1bSJed Brown 
496c4762a1bSJed Brown   PetscFunctionBegin;
4979566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.0));
498c4762a1bSJed Brown 
4999566063dSJacob Faibussowitsch   PetscCall(TSGetDM(ts, &networkdm));
5009566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localF));
5019566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localX));
5029566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localXdot));
5039566063dSJacob Faibussowitsch   PetscCall(VecSet(localF, 0.0));
504c4762a1bSJed Brown 
505c4762a1bSJed Brown   /* update ghost values of localX and localXdot */
5069566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
5079566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
508c4762a1bSJed Brown 
5099566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm, Xdot, INSERT_VALUES, localXdot));
5109566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm, Xdot, INSERT_VALUES, localXdot));
511c4762a1bSJed Brown 
5129566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localX, &xarr));
5139566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localXdot, &xdotarr));
5149566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localF, &farr));
515c4762a1bSJed Brown 
5169566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
517c4762a1bSJed Brown 
518c4762a1bSJed Brown   for (v = vStart; v < vEnd; v++) {
519c4762a1bSJed Brown     PetscInt    i, j, offsetbus, offsetgen, offsetexc, key;
520c4762a1bSJed Brown     Bus        *bus;
521c4762a1bSJed Brown     Gen        *gen;
522c4762a1bSJed Brown     Exc        *exc;
523c4762a1bSJed Brown     Load       *load;
524c4762a1bSJed Brown     PetscBool   ghostvtex;
525c4762a1bSJed Brown     PetscInt    numComps;
526c4762a1bSJed Brown     PetscScalar Yffr, Yffi; /* Real and imaginary fault admittances */
527c4762a1bSJed Brown     PetscScalar Vm, Vm2, Vm0;
528c4762a1bSJed Brown     PetscScalar Vr0 = 0, Vi0 = 0;
529c4762a1bSJed Brown     PetscScalar PD, QD;
530c4762a1bSJed Brown 
5319566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
5329566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm, v, &numComps));
533c4762a1bSJed Brown 
534c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
5359566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(networkdm, v, j, &key, &component, NULL));
536c4762a1bSJed Brown       if (key == 1) {
537c4762a1bSJed Brown         PetscInt        nconnedges;
538c4762a1bSJed Brown         const PetscInt *connedges;
539c4762a1bSJed Brown 
540c4762a1bSJed Brown         bus = (Bus *)(component);
5419566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetbus));
542c4762a1bSJed Brown         if (!ghostvtex) {
543c4762a1bSJed Brown           Vr = xarr[offsetbus];
544c4762a1bSJed Brown           Vi = xarr[offsetbus + 1];
545c4762a1bSJed Brown 
546c4762a1bSJed Brown           Yffr = bus->yff[1];
547c4762a1bSJed Brown           Yffi = bus->yff[0];
548c4762a1bSJed Brown 
549c4762a1bSJed Brown           if (user->alg_flg) {
550c4762a1bSJed Brown             Yffr += user->ybusfault[bus->id * 2 + 1];
551c4762a1bSJed Brown             Yffi += user->ybusfault[bus->id * 2];
552c4762a1bSJed Brown           }
553c4762a1bSJed Brown           Vr0 = bus->vr;
554c4762a1bSJed Brown           Vi0 = bus->vi;
555c4762a1bSJed Brown 
556c4762a1bSJed Brown           /* Network current balance residual IG + Y*V + IL = 0. Only YV is added here.
557c4762a1bSJed Brown            The generator current injection, IG, and load current injection, ID are added later
558c4762a1bSJed Brown            */
559c4762a1bSJed Brown           farr[offsetbus] += Yffi * Vr + Yffr * Vi;     /* imaginary current due to diagonal elements */
560c4762a1bSJed Brown           farr[offsetbus + 1] += Yffr * Vr - Yffi * Vi; /* real current due to diagonal elements */
561c4762a1bSJed Brown         }
562c4762a1bSJed Brown 
5639566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetSupportingEdges(networkdm, v, &nconnedges, &connedges));
564c4762a1bSJed Brown 
565c4762a1bSJed Brown         for (i = 0; i < nconnedges; i++) {
566c4762a1bSJed Brown           Branch         *branch;
567c4762a1bSJed Brown           PetscInt        keye;
568c4762a1bSJed Brown           PetscScalar     Yfti, Yftr, Vfr, Vfi, Vtr, Vti;
569c4762a1bSJed Brown           const PetscInt *cone;
570c4762a1bSJed Brown 
571c4762a1bSJed Brown           e = connedges[i];
5729566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, (void **)&branch, NULL));
573c4762a1bSJed Brown 
574c4762a1bSJed Brown           Yfti = branch->yft[0];
575c4762a1bSJed Brown           Yftr = branch->yft[1];
576c4762a1bSJed Brown 
5779566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone));
578c4762a1bSJed Brown 
579c4762a1bSJed Brown           vfrom = cone[0];
580c4762a1bSJed Brown           vto   = cone[1];
581c4762a1bSJed Brown 
5829566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, 0, &offsetfrom));
5839566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, 0, &offsetto));
584c4762a1bSJed Brown 
585c4762a1bSJed Brown           /* From bus and to bus real and imaginary voltages */
586c4762a1bSJed Brown           Vfr = xarr[offsetfrom];
587c4762a1bSJed Brown           Vfi = xarr[offsetfrom + 1];
588c4762a1bSJed Brown           Vtr = xarr[offsetto];
589c4762a1bSJed Brown           Vti = xarr[offsetto + 1];
590c4762a1bSJed Brown 
591c4762a1bSJed Brown           if (vfrom == v) {
592c4762a1bSJed Brown             farr[offsetfrom] += Yftr * Vti + Yfti * Vtr;
593c4762a1bSJed Brown             farr[offsetfrom + 1] += Yftr * Vtr - Yfti * Vti;
594c4762a1bSJed Brown           } else {
595c4762a1bSJed Brown             farr[offsetto] += Yftr * Vfi + Yfti * Vfr;
596c4762a1bSJed Brown             farr[offsetto + 1] += Yftr * Vfr - Yfti * Vfi;
597c4762a1bSJed Brown           }
598c4762a1bSJed Brown         }
599c4762a1bSJed Brown       } else if (key == 2) {
600c4762a1bSJed Brown         if (!ghostvtex) {
601c4762a1bSJed Brown           PetscScalar Eqp, Edp, delta, w; /* Generator variables */
602c4762a1bSJed Brown           PetscScalar Efd;                /* Exciter field voltage */
603c4762a1bSJed Brown           PetscScalar Id, Iq;             /* Generator dq axis currents */
604c4762a1bSJed Brown           PetscScalar IGr, IGi, Zdq_inv[4], det;
605c4762a1bSJed Brown           PetscScalar Xd, Xdp, Td0p, Xq, Xqp, Tq0p, TM, D, M, Rs; /* Generator parameters */
606c4762a1bSJed Brown 
607c4762a1bSJed Brown           gen = (Gen *)(component);
6089566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetgen));
609c4762a1bSJed Brown 
610c4762a1bSJed Brown           /* Generator state variables */
611c4762a1bSJed Brown           Eqp   = xarr[offsetgen];
612c4762a1bSJed Brown           Edp   = xarr[offsetgen + 1];
613c4762a1bSJed Brown           delta = xarr[offsetgen + 2];
614c4762a1bSJed Brown           w     = xarr[offsetgen + 3];
615c4762a1bSJed Brown           Id    = xarr[offsetgen + 4];
616c4762a1bSJed Brown           Iq    = xarr[offsetgen + 5];
617c4762a1bSJed Brown 
618c4762a1bSJed Brown           /* Generator parameters */
619c4762a1bSJed Brown           Xd   = gen->Xd;
620c4762a1bSJed Brown           Xdp  = gen->Xdp;
621c4762a1bSJed Brown           Td0p = gen->Td0p;
622c4762a1bSJed Brown           Xq   = gen->Xq;
623c4762a1bSJed Brown           Xqp  = gen->Xqp;
624c4762a1bSJed Brown           Tq0p = gen->Tq0p;
625c4762a1bSJed Brown           TM   = gen->TM;
626c4762a1bSJed Brown           D    = gen->D;
627c4762a1bSJed Brown           M    = gen->M;
628c4762a1bSJed Brown           Rs   = gen->Rs;
629c4762a1bSJed Brown 
6309566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, 2, &offsetexc));
631c4762a1bSJed Brown           Efd = xarr[offsetexc];
632c4762a1bSJed Brown 
633c4762a1bSJed Brown           /* Generator differential equations */
634c4762a1bSJed Brown           farr[offsetgen]     = (Eqp + (Xd - Xdp) * Id - Efd) / Td0p + xdotarr[offsetgen];
635c4762a1bSJed Brown           farr[offsetgen + 1] = (Edp - (Xq - Xqp) * Iq) / Tq0p + xdotarr[offsetgen + 1];
636c4762a1bSJed Brown           farr[offsetgen + 2] = -w + W_S + xdotarr[offsetgen + 2];
637c4762a1bSJed Brown           farr[offsetgen + 3] = (-TM + Edp * Id + Eqp * Iq + (Xqp - Xdp) * Id * Iq + D * (w - W_S)) / M + xdotarr[offsetgen + 3];
638c4762a1bSJed Brown 
6399566063dSJacob Faibussowitsch           PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
640c4762a1bSJed Brown 
641c4762a1bSJed Brown           /* Algebraic equations for stator currents */
642c4762a1bSJed Brown           det = Rs * Rs + Xdp * Xqp;
643c4762a1bSJed Brown 
644c4762a1bSJed Brown           Zdq_inv[0] = Rs / det;
645c4762a1bSJed Brown           Zdq_inv[1] = Xqp / det;
646c4762a1bSJed Brown           Zdq_inv[2] = -Xdp / det;
647c4762a1bSJed Brown           Zdq_inv[3] = Rs / det;
648c4762a1bSJed Brown 
649c4762a1bSJed Brown           farr[offsetgen + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
650c4762a1bSJed Brown           farr[offsetgen + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
651c4762a1bSJed Brown 
6529566063dSJacob Faibussowitsch           PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
653c4762a1bSJed Brown 
654c4762a1bSJed Brown           /* Add generator current injection to network */
655c4762a1bSJed Brown           farr[offsetbus] -= IGi;
656c4762a1bSJed Brown           farr[offsetbus + 1] -= IGr;
657c4762a1bSJed Brown         }
658c4762a1bSJed Brown       } else if (key == 3) {
659c4762a1bSJed Brown         if (!ghostvtex) {
660c4762a1bSJed Brown           PetscScalar k1, k2, KE, TE, TF, KA, KF, Vref, TA; /* Generator parameters */
661c4762a1bSJed Brown           PetscScalar Efd, RF, VR;                          /* Exciter variables */
662c4762a1bSJed Brown 
663c4762a1bSJed Brown           exc = (Exc *)(component);
6649566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetexc));
665c4762a1bSJed Brown 
666c4762a1bSJed Brown           Efd = xarr[offsetexc];
667c4762a1bSJed Brown           RF  = xarr[offsetexc + 1];
668c4762a1bSJed Brown           VR  = xarr[offsetexc + 2];
669c4762a1bSJed Brown 
670c4762a1bSJed Brown           k1   = exc->k1;
671c4762a1bSJed Brown           k2   = exc->k2;
672c4762a1bSJed Brown           KE   = exc->KE;
673c4762a1bSJed Brown           TE   = exc->TE;
674c4762a1bSJed Brown           TF   = exc->TF;
675c4762a1bSJed Brown           KA   = exc->KA;
676c4762a1bSJed Brown           KF   = exc->KF;
677c4762a1bSJed Brown           Vref = exc->Vref;
678c4762a1bSJed Brown           TA   = exc->TA;
679c4762a1bSJed Brown 
680c4762a1bSJed Brown           Vm = PetscSqrtScalar(Vd * Vd + Vq * Vq);
681c4762a1bSJed Brown           SE = k1 * PetscExpScalar(k2 * Efd);
682c4762a1bSJed Brown 
683c4762a1bSJed Brown           /* Exciter differential equations */
684c4762a1bSJed Brown           farr[offsetexc]     = (KE * Efd + SE - VR) / TE + xdotarr[offsetexc];
685c4762a1bSJed Brown           farr[offsetexc + 1] = (RF - KF * Efd / TF) / TF + xdotarr[offsetexc + 1];
686c4762a1bSJed Brown           farr[offsetexc + 2] = (VR - KA * RF + KA * KF * Efd / TF - KA * (Vref - Vm)) / TA + xdotarr[offsetexc + 2];
687c4762a1bSJed Brown         }
688c4762a1bSJed Brown       } else if (key == 4) {
689c4762a1bSJed Brown         if (!ghostvtex) {
690c4762a1bSJed Brown           PetscInt     k;
691c4762a1bSJed Brown           PetscInt     ld_nsegsp;
692c4762a1bSJed Brown           PetscInt     ld_nsegsq;
693c4762a1bSJed Brown           PetscScalar *ld_alphap;
694c4762a1bSJed Brown           PetscScalar *ld_betap, *ld_alphaq, *ld_betaq, PD0, QD0, IDr, IDi;
695c4762a1bSJed Brown 
696c4762a1bSJed Brown           load = (Load *)(component);
697c4762a1bSJed Brown 
698c4762a1bSJed Brown           /* Load Parameters */
699c4762a1bSJed Brown           ld_nsegsp = load->ld_nsegsp;
700c4762a1bSJed Brown           ld_alphap = load->ld_alphap;
701c4762a1bSJed Brown           ld_betap  = load->ld_betap;
702c4762a1bSJed Brown           ld_nsegsq = load->ld_nsegsq;
703c4762a1bSJed Brown           ld_alphaq = load->ld_alphaq;
704c4762a1bSJed Brown           ld_betaq  = load->ld_betaq;
705c4762a1bSJed Brown           PD0       = load->PD0;
706c4762a1bSJed Brown           QD0       = load->QD0;
707c4762a1bSJed Brown 
708c4762a1bSJed Brown           Vr  = xarr[offsetbus];     /* Real part of generator terminal voltage */
709c4762a1bSJed Brown           Vi  = xarr[offsetbus + 1]; /* Imaginary part of the generator terminal voltage */
710c4762a1bSJed Brown           Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
711c4762a1bSJed Brown           Vm2 = Vm * Vm;
712c4762a1bSJed Brown           Vm0 = PetscSqrtScalar(Vr0 * Vr0 + Vi0 * Vi0);
713c4762a1bSJed Brown           PD = QD = 0.0;
71457508eceSPierre Jolivet           for (k = 0; k < ld_nsegsp; k++) PD += ld_alphap[k] * PD0 * PetscPowScalar(Vm / Vm0, ld_betap[k]);
71557508eceSPierre Jolivet           for (k = 0; k < ld_nsegsq; k++) QD += ld_alphaq[k] * QD0 * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
716c4762a1bSJed Brown 
717c4762a1bSJed Brown           /* Load currents */
718c4762a1bSJed Brown           IDr = (PD * Vr + QD * Vi) / Vm2;
719c4762a1bSJed Brown           IDi = (-QD * Vr + PD * Vi) / Vm2;
720c4762a1bSJed Brown 
721c4762a1bSJed Brown           /* Load current contribution to the network */
722c4762a1bSJed Brown           farr[offsetbus] += IDi;
723c4762a1bSJed Brown           farr[offsetbus + 1] += IDr;
724c4762a1bSJed Brown         }
725c4762a1bSJed Brown       }
726c4762a1bSJed Brown     }
727c4762a1bSJed Brown   }
728c4762a1bSJed Brown 
7299566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localX, &xarr));
7309566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localXdot, &xdotarr));
7319566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localF, &farr));
7329566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localX));
7339566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localXdot));
734c4762a1bSJed Brown 
7359566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F));
7369566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F));
7379566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localF));
7383ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
739c4762a1bSJed Brown }
740c4762a1bSJed Brown 
741c4762a1bSJed Brown /* This function is used for solving the algebraic system only during fault on and
742c4762a1bSJed Brown    off times. It computes the entire F and then zeros out the part corresponding to
743c4762a1bSJed Brown    differential equations
744c4762a1bSJed Brown  F = [0;g(y)];
745c4762a1bSJed Brown */
AlgFunction(SNES snes,Vec X,Vec F,PetscCtx ctx)746*2a8381b2SBarry Smith PetscErrorCode AlgFunction(SNES snes, Vec X, Vec F, PetscCtx ctx)
747d71ae5a4SJacob Faibussowitsch {
748c4762a1bSJed Brown   DM                 networkdm;
749c4762a1bSJed Brown   Vec                localX, localF;
750c4762a1bSJed Brown   PetscInt           vfrom, vto, offsetfrom, offsetto;
751c4762a1bSJed Brown   PetscInt           v, vStart, vEnd, e;
752c4762a1bSJed Brown   PetscScalar       *farr;
753c4762a1bSJed Brown   Userctx           *user = (Userctx *)ctx;
754c4762a1bSJed Brown   const PetscScalar *xarr;
755c4762a1bSJed Brown   void              *component;
756c4762a1bSJed Brown   PetscScalar        Vr = 0, Vi = 0;
757c4762a1bSJed Brown 
758c4762a1bSJed Brown   PetscFunctionBegin;
7599566063dSJacob Faibussowitsch   PetscCall(VecSet(F, 0.0));
7609566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes, &networkdm));
7619566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localF));
7629566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm, &localX));
7639566063dSJacob Faibussowitsch   PetscCall(VecSet(localF, 0.0));
764c4762a1bSJed Brown 
765c4762a1bSJed Brown   /* update ghost values of locaX and locaXdot */
7669566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
7679566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
768c4762a1bSJed Brown 
7699566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localX, &xarr));
7709566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localF, &farr));
771c4762a1bSJed Brown 
7729566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
773c4762a1bSJed Brown 
774c4762a1bSJed Brown   for (v = vStart; v < vEnd; v++) {
775c4762a1bSJed Brown     PetscInt    i, j, offsetbus, offsetgen, key, numComps;
776c4762a1bSJed Brown     PetscScalar Yffr, Yffi, Vm, Vm2, Vm0, Vr0 = 0, Vi0 = 0, PD, QD;
777c4762a1bSJed Brown     Bus        *bus;
778c4762a1bSJed Brown     Gen        *gen;
779c4762a1bSJed Brown     Load       *load;
780c4762a1bSJed Brown     PetscBool   ghostvtex;
781c4762a1bSJed Brown 
7829566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm, v, &ghostvtex));
7839566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm, v, &numComps));
784c4762a1bSJed Brown 
785c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
7869566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(networkdm, v, j, &key, &component, NULL));
787c4762a1bSJed Brown       if (key == 1) {
788c4762a1bSJed Brown         PetscInt        nconnedges;
789c4762a1bSJed Brown         const PetscInt *connedges;
790c4762a1bSJed Brown 
791c4762a1bSJed Brown         bus = (Bus *)(component);
7929566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetbus));
793c4762a1bSJed Brown         if (!ghostvtex) {
794c4762a1bSJed Brown           Vr = xarr[offsetbus];
795c4762a1bSJed Brown           Vi = xarr[offsetbus + 1];
796c4762a1bSJed Brown 
797c4762a1bSJed Brown           Yffr = bus->yff[1];
798c4762a1bSJed Brown           Yffi = bus->yff[0];
799c4762a1bSJed Brown           if (user->alg_flg) {
800c4762a1bSJed Brown             Yffr += user->ybusfault[bus->id * 2 + 1];
801c4762a1bSJed Brown             Yffi += user->ybusfault[bus->id * 2];
802c4762a1bSJed Brown           }
803c4762a1bSJed Brown           Vr0 = bus->vr;
804c4762a1bSJed Brown           Vi0 = bus->vi;
805c4762a1bSJed Brown 
806c4762a1bSJed Brown           farr[offsetbus] += Yffi * Vr + Yffr * Vi;
807c4762a1bSJed Brown           farr[offsetbus + 1] += Yffr * Vr - Yffi * Vi;
808c4762a1bSJed Brown         }
8099566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetSupportingEdges(networkdm, v, &nconnedges, &connedges));
810c4762a1bSJed Brown 
811c4762a1bSJed Brown         for (i = 0; i < nconnedges; i++) {
812c4762a1bSJed Brown           Branch         *branch;
813c4762a1bSJed Brown           PetscInt        keye;
814c4762a1bSJed Brown           PetscScalar     Yfti, Yftr, Vfr, Vfi, Vtr, Vti;
815c4762a1bSJed Brown           const PetscInt *cone;
816c4762a1bSJed Brown 
817c4762a1bSJed Brown           e = connedges[i];
8189566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, (void **)&branch, NULL));
819c4762a1bSJed Brown 
820c4762a1bSJed Brown           Yfti = branch->yft[0];
821c4762a1bSJed Brown           Yftr = branch->yft[1];
822c4762a1bSJed Brown 
8239566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone));
824c4762a1bSJed Brown           vfrom = cone[0];
825c4762a1bSJed Brown           vto   = cone[1];
826c4762a1bSJed Brown 
8279566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, 0, &offsetfrom));
8289566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, 0, &offsetto));
829c4762a1bSJed Brown 
830c4762a1bSJed Brown           /*From bus and to bus real and imaginary voltages */
831c4762a1bSJed Brown           Vfr = xarr[offsetfrom];
832c4762a1bSJed Brown           Vfi = xarr[offsetfrom + 1];
833c4762a1bSJed Brown           Vtr = xarr[offsetto];
834c4762a1bSJed Brown           Vti = xarr[offsetto + 1];
835c4762a1bSJed Brown 
836c4762a1bSJed Brown           if (vfrom == v) {
837c4762a1bSJed Brown             farr[offsetfrom] += Yftr * Vti + Yfti * Vtr;
838c4762a1bSJed Brown             farr[offsetfrom + 1] += Yftr * Vtr - Yfti * Vti;
839c4762a1bSJed Brown           } else {
840c4762a1bSJed Brown             farr[offsetto] += Yftr * Vfi + Yfti * Vfr;
841c4762a1bSJed Brown             farr[offsetto + 1] += Yftr * Vfr - Yfti * Vfi;
842c4762a1bSJed Brown           }
843c4762a1bSJed Brown         }
844c4762a1bSJed Brown       } else if (key == 2) {
845c4762a1bSJed Brown         if (!ghostvtex) {
846c4762a1bSJed Brown           PetscScalar Eqp, Edp, delta; /* Generator variables */
847c4762a1bSJed Brown           PetscScalar Id, Iq;          /* Generator dq axis currents */
848c4762a1bSJed Brown           PetscScalar Vd, Vq, IGr, IGi, Zdq_inv[4], det;
849c4762a1bSJed Brown           PetscScalar Xdp, Xqp, Rs; /* Generator parameters */
850c4762a1bSJed Brown 
851c4762a1bSJed Brown           gen = (Gen *)(component);
8529566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetgen));
853c4762a1bSJed Brown 
854c4762a1bSJed Brown           /* Generator state variables */
855c4762a1bSJed Brown           Eqp   = xarr[offsetgen];
856c4762a1bSJed Brown           Edp   = xarr[offsetgen + 1];
857c4762a1bSJed Brown           delta = xarr[offsetgen + 2];
858c4762a1bSJed Brown           /* w     = xarr[idx+3]; not being used */
859c4762a1bSJed Brown           Id = xarr[offsetgen + 4];
860c4762a1bSJed Brown           Iq = xarr[offsetgen + 5];
861c4762a1bSJed Brown 
862c4762a1bSJed Brown           /* Generator parameters */
863c4762a1bSJed Brown           Xdp = gen->Xdp;
864c4762a1bSJed Brown           Xqp = gen->Xqp;
865c4762a1bSJed Brown           Rs  = gen->Rs;
866c4762a1bSJed Brown 
867c4762a1bSJed Brown           /* Set generator differential equation residual functions to zero */
868c4762a1bSJed Brown           farr[offsetgen]     = 0;
869c4762a1bSJed Brown           farr[offsetgen + 1] = 0;
870c4762a1bSJed Brown           farr[offsetgen + 2] = 0;
871c4762a1bSJed Brown           farr[offsetgen + 3] = 0;
872c4762a1bSJed Brown 
8739566063dSJacob Faibussowitsch           PetscCall(ri2dq(Vr, Vi, delta, &Vd, &Vq));
874c4762a1bSJed Brown 
875c4762a1bSJed Brown           /* Algebraic equations for stator currents */
876c4762a1bSJed Brown           det = Rs * Rs + Xdp * Xqp;
877c4762a1bSJed Brown 
878c4762a1bSJed Brown           Zdq_inv[0] = Rs / det;
879c4762a1bSJed Brown           Zdq_inv[1] = Xqp / det;
880c4762a1bSJed Brown           Zdq_inv[2] = -Xdp / det;
881c4762a1bSJed Brown           Zdq_inv[3] = Rs / det;
882c4762a1bSJed Brown 
883c4762a1bSJed Brown           farr[offsetgen + 4] = Zdq_inv[0] * (-Edp + Vd) + Zdq_inv[1] * (-Eqp + Vq) + Id;
884c4762a1bSJed Brown           farr[offsetgen + 5] = Zdq_inv[2] * (-Edp + Vd) + Zdq_inv[3] * (-Eqp + Vq) + Iq;
885c4762a1bSJed Brown 
886c4762a1bSJed Brown           /* Add generator current injection to network */
8879566063dSJacob Faibussowitsch           PetscCall(dq2ri(Id, Iq, delta, &IGr, &IGi));
888c4762a1bSJed Brown 
889c4762a1bSJed Brown           farr[offsetbus] -= IGi;
890c4762a1bSJed Brown           farr[offsetbus + 1] -= IGr;
891c4762a1bSJed Brown 
892c4762a1bSJed Brown           /* Vm = PetscSqrtScalar(Vd*Vd + Vq*Vq);*/ /* a compiler warning: "Value stored to 'Vm' is never read" - comment out by Hong Zhang */
893c4762a1bSJed Brown         }
894c4762a1bSJed Brown       } else if (key == 3) {
895c4762a1bSJed Brown         if (!ghostvtex) {
896c4762a1bSJed Brown           PetscInt offsetexc;
8979566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm, v, j, &offsetexc));
898c4762a1bSJed Brown           /* Set exciter differential equation residual functions equal to zero*/
899c4762a1bSJed Brown           farr[offsetexc]     = 0;
900c4762a1bSJed Brown           farr[offsetexc + 1] = 0;
901c4762a1bSJed Brown           farr[offsetexc + 2] = 0;
902c4762a1bSJed Brown         }
903c4762a1bSJed Brown       } else if (key == 4) {
904c4762a1bSJed Brown         if (!ghostvtex) {
905c4762a1bSJed Brown           PetscInt     k, ld_nsegsp, ld_nsegsq;
906c4762a1bSJed Brown           PetscScalar *ld_alphap, *ld_betap, *ld_alphaq, *ld_betaq, PD0, QD0, IDr, IDi;
907c4762a1bSJed Brown 
908c4762a1bSJed Brown           load = (Load *)(component);
909c4762a1bSJed Brown 
910c4762a1bSJed Brown           /* Load Parameters */
911c4762a1bSJed Brown           ld_nsegsp = load->ld_nsegsp;
912c4762a1bSJed Brown           ld_alphap = load->ld_alphap;
913c4762a1bSJed Brown           ld_betap  = load->ld_betap;
914c4762a1bSJed Brown           ld_nsegsq = load->ld_nsegsq;
915c4762a1bSJed Brown           ld_alphaq = load->ld_alphaq;
916c4762a1bSJed Brown           ld_betaq  = load->ld_betaq;
917c4762a1bSJed Brown 
918c4762a1bSJed Brown           PD0 = load->PD0;
919c4762a1bSJed Brown           QD0 = load->QD0;
920c4762a1bSJed Brown 
921c4762a1bSJed Brown           Vm  = PetscSqrtScalar(Vr * Vr + Vi * Vi);
922c4762a1bSJed Brown           Vm2 = Vm * Vm;
923c4762a1bSJed Brown           Vm0 = PetscSqrtScalar(Vr0 * Vr0 + Vi0 * Vi0);
924c4762a1bSJed Brown           PD = QD = 0.0;
92557508eceSPierre Jolivet           for (k = 0; k < ld_nsegsp; k++) PD += ld_alphap[k] * PD0 * PetscPowScalar(Vm / Vm0, ld_betap[k]);
92657508eceSPierre Jolivet           for (k = 0; k < ld_nsegsq; k++) QD += ld_alphaq[k] * QD0 * PetscPowScalar(Vm / Vm0, ld_betaq[k]);
927c4762a1bSJed Brown 
928c4762a1bSJed Brown           /* Load currents */
929c4762a1bSJed Brown           IDr = (PD * Vr + QD * Vi) / Vm2;
930c4762a1bSJed Brown           IDi = (-QD * Vr + PD * Vi) / Vm2;
931c4762a1bSJed Brown 
932c4762a1bSJed Brown           farr[offsetbus] += IDi;
933c4762a1bSJed Brown           farr[offsetbus + 1] += IDr;
934c4762a1bSJed Brown         }
935c4762a1bSJed Brown       }
936c4762a1bSJed Brown     }
937c4762a1bSJed Brown   }
938c4762a1bSJed Brown 
9399566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localX, &xarr));
9409566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localF, &farr));
9419566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localX));
942c4762a1bSJed Brown 
9439566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F));
9449566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F));
9459566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm, &localF));
9463ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
947c4762a1bSJed Brown }
948c4762a1bSJed Brown 
main(int argc,char ** argv)949d71ae5a4SJacob Faibussowitsch int main(int argc, char **argv)
950d71ae5a4SJacob Faibussowitsch {
951c4762a1bSJed Brown   PetscInt      i, j, *edgelist = NULL, eStart, eEnd, vStart, vEnd;
952c4762a1bSJed Brown   PetscInt      genj, excj, loadj, componentkey[5];
953c4762a1bSJed Brown   PetscInt      nc = 1; /* No. of copies (default = 1) */
954c4762a1bSJed Brown   PetscMPIInt   size, rank;
955c4762a1bSJed Brown   Vec           X, F_alg;
956c4762a1bSJed Brown   TS            ts;
957c4762a1bSJed Brown   SNES          snes_alg, snes;
958c4762a1bSJed Brown   Bus          *bus;
959c4762a1bSJed Brown   Branch       *branch;
960c4762a1bSJed Brown   Gen          *gen;
961c4762a1bSJed Brown   Exc          *exc;
962c4762a1bSJed Brown   Load         *load;
963c4762a1bSJed Brown   DM            networkdm;
964c4762a1bSJed Brown   PetscLogStage stage1;
965c4762a1bSJed Brown   Userctx       user;
966c4762a1bSJed Brown   KSP           ksp;
967c4762a1bSJed Brown   PC            pc;
968f11a936eSBarry Smith   PetscInt      numEdges = 0;
969c4762a1bSJed Brown 
970327415f7SBarry Smith   PetscFunctionBeginUser;
9719566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc, &argv, "ex9busnetworkops", help));
9729566063dSJacob Faibussowitsch   PetscCall(PetscOptionsGetInt(NULL, NULL, "-nc", &nc, NULL));
9739566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
9749566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
975c4762a1bSJed Brown 
976c4762a1bSJed Brown   /* Read initial voltage vector and Ybus */
97748a46eb9SPierre Jolivet   if (rank == 0) PetscCall(read_data(nc, &gen, &exc, &load, &bus, &branch, &edgelist));
978c4762a1bSJed Brown 
9799566063dSJacob Faibussowitsch   PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &networkdm));
9809566063dSJacob Faibussowitsch   PetscCall(DMNetworkRegisterComponent(networkdm, "branchstruct", sizeof(Branch), &componentkey[0]));
9819566063dSJacob Faibussowitsch   PetscCall(DMNetworkRegisterComponent(networkdm, "busstruct", sizeof(Bus), &componentkey[1]));
9829566063dSJacob Faibussowitsch   PetscCall(DMNetworkRegisterComponent(networkdm, "genstruct", sizeof(Gen), &componentkey[2]));
9839566063dSJacob Faibussowitsch   PetscCall(DMNetworkRegisterComponent(networkdm, "excstruct", sizeof(Exc), &componentkey[3]));
9849566063dSJacob Faibussowitsch   PetscCall(DMNetworkRegisterComponent(networkdm, "loadstruct", sizeof(Load), &componentkey[4]));
985c4762a1bSJed Brown 
9869566063dSJacob Faibussowitsch   PetscCall(PetscLogStageRegister("Create network", &stage1));
9879566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePush(stage1));
988c4762a1bSJed Brown 
989f11a936eSBarry Smith   /* Set local number of edges and edge connectivity */
990dd400576SPatrick Sanan   if (rank == 0) numEdges = NBRANCH * nc + (nc - 1);
9919566063dSJacob Faibussowitsch   PetscCall(DMNetworkSetNumSubNetworks(networkdm, PETSC_DECIDE, 1));
9929566063dSJacob Faibussowitsch   PetscCall(DMNetworkAddSubnetwork(networkdm, NULL, numEdges, edgelist, NULL));
993c4762a1bSJed Brown 
994c4762a1bSJed Brown   /* Set up the network layout */
9959566063dSJacob Faibussowitsch   PetscCall(DMNetworkLayoutSetUp(networkdm));
996c4762a1bSJed Brown 
99748a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscFree(edgelist));
998c4762a1bSJed Brown 
9992bf73ac6SHong Zhang   /* Add network components (physical parameters of nodes and branches) and number of variables */
1000dd400576SPatrick Sanan   if (rank == 0) {
10019566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetEdgeRange(networkdm, &eStart, &eEnd));
10029371c9d4SSatish Balay     genj  = 0;
10039371c9d4SSatish Balay     loadj = 0;
10049371c9d4SSatish Balay     excj  = 0;
100548a46eb9SPierre Jolivet     for (i = eStart; i < eEnd; i++) PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[0], &branch[i - eStart], 0));
1006c4762a1bSJed Brown 
10079566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetVertexRange(networkdm, &vStart, &vEnd));
1008c4762a1bSJed Brown 
1009c4762a1bSJed Brown     for (i = vStart; i < vEnd; i++) {
10109566063dSJacob Faibussowitsch       PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[1], &bus[i - vStart], 2));
1011c4762a1bSJed Brown       if (bus[i - vStart].nofgen) {
1012c4762a1bSJed Brown         for (j = 0; j < bus[i - vStart].nofgen; j++) {
1013c4762a1bSJed Brown           /* Add generator */
10149566063dSJacob Faibussowitsch           PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[2], &gen[genj++], 6));
1015c4762a1bSJed Brown           /* Add exciter */
10169566063dSJacob Faibussowitsch           PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[3], &exc[excj++], 3));
1017c4762a1bSJed Brown         }
1018c4762a1bSJed Brown       }
1019c4762a1bSJed Brown       if (bus[i - vStart].nofload) {
102048a46eb9SPierre Jolivet         for (j = 0; j < bus[i - vStart].nofload; j++) PetscCall(DMNetworkAddComponent(networkdm, i, componentkey[4], &load[loadj++], 0));
1021c4762a1bSJed Brown       }
1022c4762a1bSJed Brown     }
1023c4762a1bSJed Brown   }
1024c4762a1bSJed Brown 
10259566063dSJacob Faibussowitsch   PetscCall(DMSetUp(networkdm));
1026c4762a1bSJed Brown 
102748a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscFree5(bus, gen, load, branch, exc));
1028c4762a1bSJed Brown 
1029c4762a1bSJed Brown   /* for parallel options: Network partitioning and distribution of data */
103048a46eb9SPierre Jolivet   if (size > 1) PetscCall(DMNetworkDistribute(&networkdm, 0));
10319566063dSJacob Faibussowitsch   PetscCall(PetscLogStagePop());
1032c4762a1bSJed Brown 
10339566063dSJacob Faibussowitsch   PetscCall(DMCreateGlobalVector(networkdm, &X));
1034c4762a1bSJed Brown 
10359566063dSJacob Faibussowitsch   PetscCall(SetInitialGuess(networkdm, X));
1036c4762a1bSJed Brown 
1037c4762a1bSJed Brown   /* Options for fault simulation */
1038d0609cedSBarry Smith   PetscOptionsBegin(PETSC_COMM_WORLD, NULL, "Transient stability fault options", "");
1039c4762a1bSJed Brown   user.tfaulton  = 0.02;
1040c4762a1bSJed Brown   user.tfaultoff = 0.05;
1041c4762a1bSJed Brown   user.Rfault    = 0.0001;
1042c4762a1bSJed Brown   user.faultbus  = 8;
10439566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tfaulton", "", "", user.tfaulton, &user.tfaulton, NULL));
10449566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tfaultoff", "", "", user.tfaultoff, &user.tfaultoff, NULL));
10459566063dSJacob Faibussowitsch   PetscCall(PetscOptionsInt("-faultbus", "", "", user.faultbus, &user.faultbus, NULL));
1046c4762a1bSJed Brown   user.t0   = 0.0;
1047c4762a1bSJed Brown   user.tmax = 0.1;
10489566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-t0", "", "", user.t0, &user.t0, NULL));
10499566063dSJacob Faibussowitsch   PetscCall(PetscOptionsReal("-tmax", "", "", user.tmax, &user.tmax, NULL));
1050c4762a1bSJed Brown 
10519566063dSJacob Faibussowitsch   PetscCall(PetscMalloc1(18 * nc, &user.ybusfault));
1052ad540459SPierre Jolivet   for (i = 0; i < 18 * nc; i++) user.ybusfault[i] = 0;
1053c4762a1bSJed Brown   user.ybusfault[user.faultbus * 2 + 1] = 1 / user.Rfault;
1054d0609cedSBarry Smith   PetscOptionsEnd();
1055c4762a1bSJed Brown 
1056c4762a1bSJed Brown   /* Setup TS solver                                           */
1057c4762a1bSJed Brown   /*--------------------------------------------------------*/
10589566063dSJacob Faibussowitsch   PetscCall(TSCreate(PETSC_COMM_WORLD, &ts));
10599566063dSJacob Faibussowitsch   PetscCall(TSSetDM(ts, (DM)networkdm));
10609566063dSJacob Faibussowitsch   PetscCall(TSSetType(ts, TSCN));
1061c4762a1bSJed Brown 
10629566063dSJacob Faibussowitsch   PetscCall(TSGetSNES(ts, &snes));
10639566063dSJacob Faibussowitsch   PetscCall(SNESGetKSP(snes, &ksp));
10649566063dSJacob Faibussowitsch   PetscCall(KSPGetPC(ksp, &pc));
10659566063dSJacob Faibussowitsch   PetscCall(PCSetType(pc, PCBJACOBI));
1066c4762a1bSJed Brown 
10678434afd1SBarry Smith   PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)FormIFunction, &user));
10689566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, user.tfaulton));
10699566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
10709566063dSJacob Faibussowitsch   PetscCall(TSSetTimeStep(ts, 0.01));
10719566063dSJacob Faibussowitsch   PetscCall(TSSetFromOptions(ts));
1072c4762a1bSJed Brown 
1073c4762a1bSJed Brown   /*user.alg_flg = PETSC_TRUE is the period when fault exists. We add fault admittance to Ybus matrix.
1074c4762a1bSJed Brown     eg, fault bus is 8. Y88(new)=Y88(old)+Yfault. */
1075c4762a1bSJed Brown   user.alg_flg = PETSC_FALSE;
1076c4762a1bSJed Brown 
1077c4762a1bSJed Brown   /* Prefault period */
107848a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "... (1) Prefault period ... \n"));
1079c4762a1bSJed Brown 
10809566063dSJacob Faibussowitsch   PetscCall(TSSetSolution(ts, X));
10819566063dSJacob Faibussowitsch   PetscCall(TSSetUp(ts));
10829566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
1083c4762a1bSJed Brown 
1084c4762a1bSJed Brown   /* Create the nonlinear solver for solving the algebraic system */
10859566063dSJacob Faibussowitsch   PetscCall(VecDuplicate(X, &F_alg));
1086c4762a1bSJed Brown 
10879566063dSJacob Faibussowitsch   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_alg));
10889566063dSJacob Faibussowitsch   PetscCall(SNESSetDM(snes_alg, (DM)networkdm));
10899566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, &user));
10909566063dSJacob Faibussowitsch   PetscCall(SNESSetOptionsPrefix(snes_alg, "alg_"));
10919566063dSJacob Faibussowitsch   PetscCall(SNESSetFromOptions(snes_alg));
1092c4762a1bSJed Brown 
1093c4762a1bSJed Brown   /* Apply disturbance - resistive fault at user.faultbus */
1094c4762a1bSJed Brown   /* This is done by adding shunt conductance to the diagonal location
1095c4762a1bSJed Brown      in the Ybus matrix */
1096c4762a1bSJed Brown   user.alg_flg = PETSC_TRUE;
1097c4762a1bSJed Brown 
1098c4762a1bSJed Brown   /* Solve the algebraic equations */
109948a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n... (2) Apply disturbance, solve algebraic equations ... \n"));
11009566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes_alg, NULL, X));
1101c4762a1bSJed Brown 
1102c4762a1bSJed Brown   /* Disturbance period */
11039566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, user.tfaulton));
11049566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, user.tfaultoff));
11059566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
11068434afd1SBarry Smith   PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)FormIFunction, &user));
1107c4762a1bSJed Brown 
1108c4762a1bSJed Brown   user.alg_flg = PETSC_TRUE;
110948a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n... (3) Disturbance period ... \n"));
11109566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
1111c4762a1bSJed Brown 
1112c4762a1bSJed Brown   /* Remove the fault */
11139566063dSJacob Faibussowitsch   PetscCall(SNESSetFunction(snes_alg, F_alg, AlgFunction, &user));
1114c4762a1bSJed Brown 
1115c4762a1bSJed Brown   user.alg_flg = PETSC_FALSE;
1116c4762a1bSJed Brown   /* Solve the algebraic equations */
111748a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n... (4) Remove fault, solve algebraic equations ... \n"));
11189566063dSJacob Faibussowitsch   PetscCall(SNESSolve(snes_alg, NULL, X));
11199566063dSJacob Faibussowitsch   PetscCall(SNESDestroy(&snes_alg));
1120c4762a1bSJed Brown 
1121c4762a1bSJed Brown   /* Post-disturbance period */
11229566063dSJacob Faibussowitsch   PetscCall(TSSetTime(ts, user.tfaultoff));
11239566063dSJacob Faibussowitsch   PetscCall(TSSetMaxTime(ts, user.tmax));
11249566063dSJacob Faibussowitsch   PetscCall(TSSetExactFinalTime(ts, TS_EXACTFINALTIME_STEPOVER));
11258434afd1SBarry Smith   PetscCall(TSSetIFunction(ts, NULL, (TSIFunctionFn *)FormIFunction, &user));
1126c4762a1bSJed Brown 
1127c4762a1bSJed Brown   user.alg_flg = PETSC_FALSE;
112848a46eb9SPierre Jolivet   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "\n... (5) Post-disturbance period ... \n"));
11299566063dSJacob Faibussowitsch   PetscCall(TSSolve(ts, X));
1130c4762a1bSJed Brown 
11319566063dSJacob Faibussowitsch   PetscCall(PetscFree(user.ybusfault));
11329566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&F_alg));
11339566063dSJacob Faibussowitsch   PetscCall(VecDestroy(&X));
11349566063dSJacob Faibussowitsch   PetscCall(DMDestroy(&networkdm));
11359566063dSJacob Faibussowitsch   PetscCall(TSDestroy(&ts));
11369566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
1137b122ec5aSJacob Faibussowitsch   return 0;
1138c4762a1bSJed Brown }
1139c4762a1bSJed Brown 
1140c4762a1bSJed Brown /*TEST
1141c4762a1bSJed Brown 
1142c4762a1bSJed Brown    build:
1143dfd57a17SPierre Jolivet       requires: double !complex !defined(PETSC_USE_64BIT_INDICES)
1144c4762a1bSJed Brown 
1145c4762a1bSJed Brown    test:
1146c4762a1bSJed Brown       args: -ts_monitor -snes_converged_reason -alg_snes_converged_reason
1147c4762a1bSJed Brown       localrunfiles: X.bin Ybus.bin ex9busnetworkops
1148c4762a1bSJed Brown 
1149c4762a1bSJed Brown    test:
1150c4762a1bSJed Brown       suffix: 2
1151c4762a1bSJed Brown       nsize: 2
1152c4762a1bSJed Brown       args: -ts_monitor -snes_converged_reason -alg_snes_converged_reason
1153c4762a1bSJed Brown       localrunfiles: X.bin Ybus.bin ex9busnetworkops
1154c4762a1bSJed Brown 
1155c4762a1bSJed Brown TEST*/
1156