1c4762a1bSJed Brown /* function subroutines used by power.c */
2c4762a1bSJed Brown
3c4762a1bSJed Brown #include "power.h"
4c4762a1bSJed Brown
GetListofEdges_Power(PFDATA * pfdata,PetscInt * edgelist)5d71ae5a4SJacob Faibussowitsch PetscErrorCode GetListofEdges_Power(PFDATA *pfdata, PetscInt *edgelist)
6d71ae5a4SJacob Faibussowitsch {
7c4762a1bSJed Brown PetscInt i, fbus, tbus, nbranches = pfdata->nbranch;
8c4762a1bSJed Brown EDGE_Power branch = pfdata->branch;
9c4762a1bSJed Brown PetscBool netview = PETSC_FALSE;
10c4762a1bSJed Brown
11c4762a1bSJed Brown PetscFunctionBegin;
129566063dSJacob Faibussowitsch PetscCall(PetscOptionsHasName(NULL, NULL, "-powernet_view", &netview));
13c4762a1bSJed Brown for (i = 0; i < nbranches; i++) {
14c4762a1bSJed Brown fbus = branch[i].internal_i;
15c4762a1bSJed Brown tbus = branch[i].internal_j;
16c4762a1bSJed Brown edgelist[2 * i] = fbus;
17c4762a1bSJed Brown edgelist[2 * i + 1] = tbus;
1848a46eb9SPierre Jolivet if (netview) PetscCall(PetscPrintf(PETSC_COMM_SELF, "branch %" PetscInt_FMT ", bus[%" PetscInt_FMT "] -> bus[%" PetscInt_FMT "]\n", i, fbus, tbus));
19c4762a1bSJed Brown }
20c4762a1bSJed Brown if (netview) {
21c4762a1bSJed Brown for (i = 0; i < pfdata->nbus; i++) {
22c4762a1bSJed Brown if (pfdata->bus[i].ngen) {
2363a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " bus %" PetscInt_FMT ": gen\n", i));
24c4762a1bSJed Brown } else if (pfdata->bus[i].nload) {
2563a3b9bcSJacob Faibussowitsch PetscCall(PetscPrintf(PETSC_COMM_SELF, " bus %" PetscInt_FMT ": load\n", i));
26c4762a1bSJed Brown }
27c4762a1bSJed Brown }
28c4762a1bSJed Brown }
293ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
30c4762a1bSJed Brown }
31c4762a1bSJed Brown
FormJacobian_Power_private(DM networkdm,Vec localX,Mat J,PetscInt nv,PetscInt ne,const PetscInt * vtx,const PetscInt * edges,void * appctx)32d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian_Power_private(DM networkdm, Vec localX, Mat J, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
33d71ae5a4SJacob Faibussowitsch {
34c4762a1bSJed Brown const PetscScalar *xarr;
35c4762a1bSJed Brown PetscInt i, v, row[2], col[8], e, vfrom, vto;
36c4762a1bSJed Brown PetscInt offsetfrom, offsetto, goffsetfrom, goffsetto, numComps;
37c4762a1bSJed Brown PetscScalar values[8];
38c4762a1bSJed Brown PetscInt j, key, offset, goffset;
39c4762a1bSJed Brown PetscScalar Vm;
40c4762a1bSJed Brown UserCtx_Power *user_power = (UserCtx_Power *)appctx;
41c4762a1bSJed Brown PetscScalar Sbase = user_power->Sbase;
42c4762a1bSJed Brown VERTEX_Power bus;
43c4762a1bSJed Brown PetscBool ghostvtex;
44c4762a1bSJed Brown void *component;
45c4762a1bSJed Brown
46c4762a1bSJed Brown PetscFunctionBegin;
479566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr));
48c4762a1bSJed Brown
49c4762a1bSJed Brown for (v = 0; v < nv; v++) {
509566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghostvtex));
51c4762a1bSJed Brown
529566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &numComps));
53c4762a1bSJed Brown for (j = 0; j < numComps; j++) {
549566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &offset));
559566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &goffset));
569566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[v], j, &key, &component, NULL));
57c4762a1bSJed Brown
58c4762a1bSJed Brown if (key == user_power->compkey_bus) {
59c4762a1bSJed Brown PetscInt nconnedges;
60c4762a1bSJed Brown const PetscInt *connedges;
61c4762a1bSJed Brown
62*57508eceSPierre Jolivet bus = (VERTEX_Power)component;
63c4762a1bSJed Brown if (!ghostvtex) {
64c4762a1bSJed Brown /* Handle reference bus constrained dofs */
65c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
669371c9d4SSatish Balay row[0] = goffset;
679371c9d4SSatish Balay row[1] = goffset + 1;
689371c9d4SSatish Balay col[0] = goffset;
699371c9d4SSatish Balay col[1] = goffset + 1;
709371c9d4SSatish Balay col[2] = goffset;
719371c9d4SSatish Balay col[3] = goffset + 1;
729371c9d4SSatish Balay values[0] = 1.0;
739371c9d4SSatish Balay values[1] = 0.0;
749371c9d4SSatish Balay values[2] = 0.0;
759371c9d4SSatish Balay values[3] = 1.0;
769566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, row, 2, col, values, ADD_VALUES));
77c4762a1bSJed Brown break;
78c4762a1bSJed Brown }
79c4762a1bSJed Brown
80c4762a1bSJed Brown Vm = xarr[offset + 1];
81c4762a1bSJed Brown
82c4762a1bSJed Brown /* Shunt injections */
839371c9d4SSatish Balay row[0] = goffset;
849371c9d4SSatish Balay row[1] = goffset + 1;
859371c9d4SSatish Balay col[0] = goffset;
869371c9d4SSatish Balay col[1] = goffset + 1;
87c4762a1bSJed Brown values[0] = values[1] = values[2] = values[3] = 0.0;
88c4762a1bSJed Brown if (bus->ide != PV_BUS) {
89c4762a1bSJed Brown values[1] = 2.0 * Vm * bus->gl / Sbase;
90c4762a1bSJed Brown values[3] = -2.0 * Vm * bus->bl / Sbase;
91c4762a1bSJed Brown }
929566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 2, row, 2, col, values, ADD_VALUES));
93c4762a1bSJed Brown }
94c4762a1bSJed Brown
959566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges));
96c4762a1bSJed Brown
97c4762a1bSJed Brown for (i = 0; i < nconnedges; i++) {
98c4762a1bSJed Brown EDGE_Power branch;
99c4762a1bSJed Brown VERTEX_Power busf, bust;
100c4762a1bSJed Brown PetscInt keyf, keyt;
101c4762a1bSJed Brown PetscScalar Gff, Bff, Gft, Bft, Gtf, Btf, Gtt, Btt;
102c4762a1bSJed Brown const PetscInt *cone;
103c4762a1bSJed Brown PetscScalar Vmf, Vmt, thetaf, thetat, thetaft, thetatf;
104c4762a1bSJed Brown
105c4762a1bSJed Brown e = connedges[i];
1069566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &key, (void **)&branch, NULL));
107c4762a1bSJed Brown if (!branch->status) continue;
108c4762a1bSJed Brown
109c4762a1bSJed Brown Gff = branch->yff[0];
110c4762a1bSJed Brown Bff = branch->yff[1];
111c4762a1bSJed Brown Gft = branch->yft[0];
112c4762a1bSJed Brown Bft = branch->yft[1];
113c4762a1bSJed Brown Gtf = branch->ytf[0];
114c4762a1bSJed Brown Btf = branch->ytf[1];
115c4762a1bSJed Brown Gtt = branch->ytt[0];
116c4762a1bSJed Brown Btt = branch->ytt[1];
117c4762a1bSJed Brown
1189566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, edges[e], &cone));
119c4762a1bSJed Brown vfrom = cone[0];
120c4762a1bSJed Brown vto = cone[1];
121c4762a1bSJed Brown
1229566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom));
1239566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto));
1249566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &goffsetfrom));
1259566063dSJacob Faibussowitsch PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vto, ALL_COMPONENTS, &goffsetto));
126c4762a1bSJed Brown
127c4762a1bSJed Brown if (goffsetto < 0) goffsetto = -goffsetto - 1;
128c4762a1bSJed Brown
129c4762a1bSJed Brown thetaf = xarr[offsetfrom];
130c4762a1bSJed Brown Vmf = xarr[offsetfrom + 1];
131c4762a1bSJed Brown thetat = xarr[offsetto];
132c4762a1bSJed Brown Vmt = xarr[offsetto + 1];
133c4762a1bSJed Brown thetaft = thetaf - thetat;
134c4762a1bSJed Brown thetatf = thetat - thetaf;
135c4762a1bSJed Brown
1369566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vfrom, 0, &keyf, (void **)&busf, NULL));
1379566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vto, 0, &keyt, (void **)&bust, NULL));
138c4762a1bSJed Brown
139c4762a1bSJed Brown if (vfrom == vtx[v]) {
140c4762a1bSJed Brown if (busf->ide != REF_BUS) {
141c4762a1bSJed Brown /* farr[offsetfrom] += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); */
142c4762a1bSJed Brown row[0] = goffsetfrom;
1439371c9d4SSatish Balay col[0] = goffsetfrom;
1449371c9d4SSatish Balay col[1] = goffsetfrom + 1;
1459371c9d4SSatish Balay col[2] = goffsetto;
1469371c9d4SSatish Balay col[3] = goffsetto + 1;
147c4762a1bSJed Brown values[0] = Vmf * Vmt * (Gft * -PetscSinScalar(thetaft) + Bft * PetscCosScalar(thetaft)); /* df_dthetaf */
148c4762a1bSJed Brown values[1] = 2.0 * Gff * Vmf + Vmt * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); /* df_dVmf */
149c4762a1bSJed Brown values[2] = Vmf * Vmt * (Gft * PetscSinScalar(thetaft) + Bft * -PetscCosScalar(thetaft)); /* df_dthetat */
150c4762a1bSJed Brown values[3] = Vmf * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft)); /* df_dVmt */
151c4762a1bSJed Brown
1529566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES));
153c4762a1bSJed Brown }
154c4762a1bSJed Brown if (busf->ide != PV_BUS && busf->ide != REF_BUS) {
155c4762a1bSJed Brown row[0] = goffsetfrom + 1;
1569371c9d4SSatish Balay col[0] = goffsetfrom;
1579371c9d4SSatish Balay col[1] = goffsetfrom + 1;
1589371c9d4SSatish Balay col[2] = goffsetto;
1599371c9d4SSatish Balay col[3] = goffsetto + 1;
160c4762a1bSJed Brown /* farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */
161c4762a1bSJed Brown values[0] = Vmf * Vmt * (Bft * PetscSinScalar(thetaft) + Gft * PetscCosScalar(thetaft));
162c4762a1bSJed Brown values[1] = -2.0 * Bff * Vmf + Vmt * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft));
163c4762a1bSJed Brown values[2] = Vmf * Vmt * (-Bft * PetscSinScalar(thetaft) + Gft * -PetscCosScalar(thetaft));
164c4762a1bSJed Brown values[3] = Vmf * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft));
165c4762a1bSJed Brown
1669566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES));
167c4762a1bSJed Brown }
168c4762a1bSJed Brown } else {
169c4762a1bSJed Brown if (bust->ide != REF_BUS) {
170c4762a1bSJed Brown row[0] = goffsetto;
1719371c9d4SSatish Balay col[0] = goffsetto;
1729371c9d4SSatish Balay col[1] = goffsetto + 1;
1739371c9d4SSatish Balay col[2] = goffsetfrom;
1749371c9d4SSatish Balay col[3] = goffsetfrom + 1;
175c4762a1bSJed Brown /* farr[offsetto] += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */
176c4762a1bSJed Brown values[0] = Vmt * Vmf * (Gtf * -PetscSinScalar(thetatf) + Btf * PetscCosScalar(thetaft)); /* df_dthetat */
177c4762a1bSJed Brown values[1] = 2.0 * Gtt * Vmt + Vmf * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); /* df_dVmt */
178c4762a1bSJed Brown values[2] = Vmt * Vmf * (Gtf * PetscSinScalar(thetatf) + Btf * -PetscCosScalar(thetatf)); /* df_dthetaf */
179c4762a1bSJed Brown values[3] = Vmt * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf)); /* df_dVmf */
180c4762a1bSJed Brown
1819566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES));
182c4762a1bSJed Brown }
183c4762a1bSJed Brown if (bust->ide != PV_BUS && bust->ide != REF_BUS) {
184c4762a1bSJed Brown row[0] = goffsetto + 1;
1859371c9d4SSatish Balay col[0] = goffsetto;
1869371c9d4SSatish Balay col[1] = goffsetto + 1;
1879371c9d4SSatish Balay col[2] = goffsetfrom;
1889371c9d4SSatish Balay col[3] = goffsetfrom + 1;
189c4762a1bSJed Brown /* farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */
190c4762a1bSJed Brown values[0] = Vmt * Vmf * (Btf * PetscSinScalar(thetatf) + Gtf * PetscCosScalar(thetatf));
191c4762a1bSJed Brown values[1] = -2.0 * Btt * Vmt + Vmf * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf));
192c4762a1bSJed Brown values[2] = Vmt * Vmf * (-Btf * PetscSinScalar(thetatf) + Gtf * -PetscCosScalar(thetatf));
193c4762a1bSJed Brown values[3] = Vmt * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf));
194c4762a1bSJed Brown
1959566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 4, col, values, ADD_VALUES));
196c4762a1bSJed Brown }
197c4762a1bSJed Brown }
198c4762a1bSJed Brown }
199c4762a1bSJed Brown if (!ghostvtex && bus->ide == PV_BUS) {
2009371c9d4SSatish Balay row[0] = goffset + 1;
2019371c9d4SSatish Balay col[0] = goffset + 1;
2029371c9d4SSatish Balay values[0] = 1.0;
203c4762a1bSJed Brown if (user_power->jac_error) values[0] = 50.0;
2049566063dSJacob Faibussowitsch PetscCall(MatSetValues(J, 1, row, 1, col, values, ADD_VALUES));
205c4762a1bSJed Brown }
206c4762a1bSJed Brown }
207c4762a1bSJed Brown }
208c4762a1bSJed Brown }
209c4762a1bSJed Brown
2109566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr));
2113ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
212c4762a1bSJed Brown }
213c4762a1bSJed Brown
FormJacobian_Power(SNES snes,Vec X,Mat J,Mat Jpre,void * appctx)214d71ae5a4SJacob Faibussowitsch PetscErrorCode FormJacobian_Power(SNES snes, Vec X, Mat J, Mat Jpre, void *appctx)
215d71ae5a4SJacob Faibussowitsch {
216c4762a1bSJed Brown DM networkdm;
217c4762a1bSJed Brown Vec localX;
218c4762a1bSJed Brown PetscInt nv, ne;
219c4762a1bSJed Brown const PetscInt *vtx, *edges;
220c4762a1bSJed Brown
221c4762a1bSJed Brown PetscFunctionBegin;
2229566063dSJacob Faibussowitsch PetscCall(MatZeroEntries(J));
223c4762a1bSJed Brown
2249566063dSJacob Faibussowitsch PetscCall(SNESGetDM(snes, &networkdm));
2259566063dSJacob Faibussowitsch PetscCall(DMGetLocalVector(networkdm, &localX));
226c4762a1bSJed Brown
2279566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
2289566063dSJacob Faibussowitsch PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
229c4762a1bSJed Brown
2309566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
2319566063dSJacob Faibussowitsch PetscCall(FormJacobian_Power_private(networkdm, localX, J, nv, ne, vtx, edges, appctx));
232c4762a1bSJed Brown
2339566063dSJacob Faibussowitsch PetscCall(DMRestoreLocalVector(networkdm, &localX));
234c4762a1bSJed Brown
2359566063dSJacob Faibussowitsch PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
2369566063dSJacob Faibussowitsch PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
2373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
238c4762a1bSJed Brown }
239c4762a1bSJed Brown
FormFunction_Power(DM networkdm,Vec localX,Vec localF,PetscInt nv,PetscInt ne,const PetscInt * vtx,const PetscInt * edges,void * appctx)240d71ae5a4SJacob Faibussowitsch PetscErrorCode FormFunction_Power(DM networkdm, Vec localX, Vec localF, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
241d71ae5a4SJacob Faibussowitsch {
242c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power *)appctx;
243c4762a1bSJed Brown PetscInt e, v, vfrom, vto;
244c4762a1bSJed Brown const PetscScalar *xarr;
245c4762a1bSJed Brown PetscScalar *farr;
246c4762a1bSJed Brown PetscInt offsetfrom, offsetto, offset, i, j, key, numComps;
247c4762a1bSJed Brown PetscScalar Vm;
248c4762a1bSJed Brown PetscScalar Sbase = User->Sbase;
249c4762a1bSJed Brown VERTEX_Power bus = NULL;
250c4762a1bSJed Brown GEN gen;
251c4762a1bSJed Brown LOAD load;
252c4762a1bSJed Brown PetscBool ghostvtex;
253c4762a1bSJed Brown void *component;
254c4762a1bSJed Brown
255c4762a1bSJed Brown PetscFunctionBegin;
2569566063dSJacob Faibussowitsch PetscCall(VecGetArrayRead(localX, &xarr));
2579566063dSJacob Faibussowitsch PetscCall(VecGetArray(localF, &farr));
258c4762a1bSJed Brown
259c4762a1bSJed Brown for (v = 0; v < nv; v++) {
2609566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghostvtex));
2619566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &numComps));
2629566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[v], ALL_COMPONENTS, &offset));
263c4762a1bSJed Brown
264c4762a1bSJed Brown for (j = 0; j < numComps; j++) {
2659566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[v], j, &key, &component, NULL));
266c4762a1bSJed Brown if (key == User->compkey_bus) {
267c4762a1bSJed Brown PetscInt nconnedges;
268c4762a1bSJed Brown const PetscInt *connedges;
269c4762a1bSJed Brown
270*57508eceSPierre Jolivet bus = (VERTEX_Power)component;
271c4762a1bSJed Brown /* Handle reference bus constrained dofs */
272c4762a1bSJed Brown if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
273c4762a1bSJed Brown farr[offset] = xarr[offset] - bus->va * PETSC_PI / 180.0;
274c4762a1bSJed Brown farr[offset + 1] = xarr[offset + 1] - bus->vm;
275c4762a1bSJed Brown break;
276c4762a1bSJed Brown }
277c4762a1bSJed Brown
278c4762a1bSJed Brown if (!ghostvtex) {
279c4762a1bSJed Brown Vm = xarr[offset + 1];
280c4762a1bSJed Brown
281c4762a1bSJed Brown /* Shunt injections */
282c4762a1bSJed Brown farr[offset] += Vm * Vm * bus->gl / Sbase;
283c4762a1bSJed Brown if (bus->ide != PV_BUS) farr[offset + 1] += -Vm * Vm * bus->bl / Sbase;
284c4762a1bSJed Brown }
285c4762a1bSJed Brown
2869566063dSJacob Faibussowitsch PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges));
287c4762a1bSJed Brown for (i = 0; i < nconnedges; i++) {
288c4762a1bSJed Brown EDGE_Power branch;
289c4762a1bSJed Brown PetscInt keye;
290c4762a1bSJed Brown PetscScalar Gff, Bff, Gft, Bft, Gtf, Btf, Gtt, Btt;
291c4762a1bSJed Brown const PetscInt *cone;
292c4762a1bSJed Brown PetscScalar Vmf, Vmt, thetaf, thetat, thetaft, thetatf;
293c4762a1bSJed Brown
294c4762a1bSJed Brown e = connedges[i];
2959566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, (void **)&branch, NULL));
296c4762a1bSJed Brown if (!branch->status) continue;
297c4762a1bSJed Brown Gff = branch->yff[0];
298c4762a1bSJed Brown Bff = branch->yff[1];
299c4762a1bSJed Brown Gft = branch->yft[0];
300c4762a1bSJed Brown Bft = branch->yft[1];
301c4762a1bSJed Brown Gtf = branch->ytf[0];
302c4762a1bSJed Brown Btf = branch->ytf[1];
303c4762a1bSJed Brown Gtt = branch->ytt[0];
304c4762a1bSJed Brown Btt = branch->ytt[1];
305c4762a1bSJed Brown
3069566063dSJacob Faibussowitsch PetscCall(DMNetworkGetConnectedVertices(networkdm, e, &cone));
307c4762a1bSJed Brown vfrom = cone[0];
308c4762a1bSJed Brown vto = cone[1];
309c4762a1bSJed Brown
3109566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vfrom, ALL_COMPONENTS, &offsetfrom));
3119566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vto, ALL_COMPONENTS, &offsetto));
312c4762a1bSJed Brown
313c4762a1bSJed Brown thetaf = xarr[offsetfrom];
314c4762a1bSJed Brown Vmf = xarr[offsetfrom + 1];
315c4762a1bSJed Brown thetat = xarr[offsetto];
316c4762a1bSJed Brown Vmt = xarr[offsetto + 1];
317c4762a1bSJed Brown thetaft = thetaf - thetat;
318c4762a1bSJed Brown thetatf = thetat - thetaf;
319c4762a1bSJed Brown
320c4762a1bSJed Brown if (vfrom == vtx[v]) {
321c4762a1bSJed Brown farr[offsetfrom] += Gff * Vmf * Vmf + Vmf * Vmt * (Gft * PetscCosScalar(thetaft) + Bft * PetscSinScalar(thetaft));
322c4762a1bSJed Brown farr[offsetfrom + 1] += -Bff * Vmf * Vmf + Vmf * Vmt * (-Bft * PetscCosScalar(thetaft) + Gft * PetscSinScalar(thetaft));
323c4762a1bSJed Brown } else {
324c4762a1bSJed Brown farr[offsetto] += Gtt * Vmt * Vmt + Vmt * Vmf * (Gtf * PetscCosScalar(thetatf) + Btf * PetscSinScalar(thetatf));
325c4762a1bSJed Brown farr[offsetto + 1] += -Btt * Vmt * Vmt + Vmt * Vmf * (-Btf * PetscCosScalar(thetatf) + Gtf * PetscSinScalar(thetatf));
326c4762a1bSJed Brown }
327c4762a1bSJed Brown }
328c4762a1bSJed Brown } else if (key == User->compkey_gen) {
329c4762a1bSJed Brown if (!ghostvtex) {
330*57508eceSPierre Jolivet gen = (GEN)component;
331c4762a1bSJed Brown if (!gen->status) continue;
332c4762a1bSJed Brown farr[offset] += -gen->pg / Sbase;
333c4762a1bSJed Brown farr[offset + 1] += -gen->qg / Sbase;
334c4762a1bSJed Brown }
335c4762a1bSJed Brown } else if (key == User->compkey_load) {
336c4762a1bSJed Brown if (!ghostvtex) {
337*57508eceSPierre Jolivet load = (LOAD)component;
338c4762a1bSJed Brown farr[offset] += load->pl / Sbase;
339c4762a1bSJed Brown farr[offset + 1] += load->ql / Sbase;
340c4762a1bSJed Brown }
341c4762a1bSJed Brown }
342c4762a1bSJed Brown }
343ad540459SPierre Jolivet if (bus && bus->ide == PV_BUS) farr[offset + 1] = xarr[offset + 1] - bus->vm;
344c4762a1bSJed Brown }
3459566063dSJacob Faibussowitsch PetscCall(VecRestoreArrayRead(localX, &xarr));
3469566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localF, &farr));
3473ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
348c4762a1bSJed Brown }
349c4762a1bSJed Brown
SetInitialGuess_Power(DM networkdm,Vec localX,PetscInt nv,PetscInt ne,const PetscInt * vtx,const PetscInt * edges,void * appctx)350d71ae5a4SJacob Faibussowitsch PetscErrorCode SetInitialGuess_Power(DM networkdm, Vec localX, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
351d71ae5a4SJacob Faibussowitsch {
352c4762a1bSJed Brown VERTEX_Power bus;
353c4762a1bSJed Brown PetscInt i;
354c4762a1bSJed Brown GEN gen;
3552bf73ac6SHong Zhang PetscBool ghostvtex, sharedv;
356c4762a1bSJed Brown PetscScalar *xarr;
357c4762a1bSJed Brown PetscInt key, numComps, j, offset;
358c4762a1bSJed Brown void *component;
359c4762a1bSJed Brown PetscMPIInt rank;
360c4762a1bSJed Brown MPI_Comm comm;
361c4762a1bSJed Brown UserCtx_Power *User = (UserCtx_Power *)appctx;
362c4762a1bSJed Brown
363c4762a1bSJed Brown PetscFunctionBegin;
3649566063dSJacob Faibussowitsch PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm));
3659566063dSJacob Faibussowitsch PetscCallMPI(MPI_Comm_rank(comm, &rank));
3669566063dSJacob Faibussowitsch PetscCall(VecGetArray(localX, &xarr));
367c4762a1bSJed Brown for (i = 0; i < nv; i++) {
3689566063dSJacob Faibussowitsch PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex));
3699566063dSJacob Faibussowitsch PetscCall(DMNetworkIsSharedVertex(networkdm, vtx[i], &sharedv));
3702bf73ac6SHong Zhang if (ghostvtex || sharedv) continue;
371c4762a1bSJed Brown
3729566063dSJacob Faibussowitsch PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset));
3739566063dSJacob Faibussowitsch PetscCall(DMNetworkGetNumComponents(networkdm, vtx[i], &numComps));
374c4762a1bSJed Brown for (j = 0; j < numComps; j++) {
3759566063dSJacob Faibussowitsch PetscCall(DMNetworkGetComponent(networkdm, vtx[i], j, &key, &component, NULL));
376c4762a1bSJed Brown if (key == User->compkey_bus) {
377*57508eceSPierre Jolivet bus = (VERTEX_Power)component;
378c4762a1bSJed Brown xarr[offset] = bus->va * PETSC_PI / 180.0;
379c4762a1bSJed Brown xarr[offset + 1] = bus->vm;
380c4762a1bSJed Brown } else if (key == User->compkey_gen) {
381*57508eceSPierre Jolivet gen = (GEN)component;
382c4762a1bSJed Brown if (!gen->status) continue;
383c4762a1bSJed Brown xarr[offset + 1] = gen->vs;
384c4762a1bSJed Brown break;
385c4762a1bSJed Brown }
386c4762a1bSJed Brown }
387c4762a1bSJed Brown }
3889566063dSJacob Faibussowitsch PetscCall(VecRestoreArray(localX, &xarr));
3893ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS);
390c4762a1bSJed Brown }
391