xref: /petsc/src/snes/tutorials/network/power/pffunctions.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown /* function subroutines used by power.c */
2*c4762a1bSJed Brown 
3*c4762a1bSJed Brown #include "power.h"
4*c4762a1bSJed Brown 
5*c4762a1bSJed Brown PetscErrorCode GetListofEdges_Power(PFDATA *pfdata,PetscInt *edgelist)
6*c4762a1bSJed Brown {
7*c4762a1bSJed Brown   PetscErrorCode ierr;
8*c4762a1bSJed Brown   PetscInt       i,fbus,tbus,nbranches=pfdata->nbranch;
9*c4762a1bSJed Brown   EDGE_Power     branch=pfdata->branch;
10*c4762a1bSJed Brown   PetscBool      netview=PETSC_FALSE;
11*c4762a1bSJed Brown 
12*c4762a1bSJed Brown   PetscFunctionBegin;
13*c4762a1bSJed Brown   ierr = PetscOptionsHasName(NULL,NULL, "-powernet_view",&netview);CHKERRQ(ierr);
14*c4762a1bSJed Brown   for (i=0; i<nbranches; i++) {
15*c4762a1bSJed Brown     fbus = branch[i].internal_i;
16*c4762a1bSJed Brown     tbus = branch[i].internal_j;
17*c4762a1bSJed Brown     edgelist[2*i]   = fbus;
18*c4762a1bSJed Brown     edgelist[2*i+1] = tbus;
19*c4762a1bSJed Brown     if (netview) {
20*c4762a1bSJed Brown       ierr = PetscPrintf(PETSC_COMM_SELF,"branch %d, bus[%d] -> bus[%d]\n",i,fbus,tbus);CHKERRQ(ierr);
21*c4762a1bSJed Brown     }
22*c4762a1bSJed Brown   }
23*c4762a1bSJed Brown   if (netview) {
24*c4762a1bSJed Brown     for (i=0; i<pfdata->nbus; i++) {
25*c4762a1bSJed Brown       if (pfdata->bus[i].ngen) {
26*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF," bus %D: gen\n",i);CHKERRQ(ierr);
27*c4762a1bSJed Brown       } else if (pfdata->bus[i].nload) {
28*c4762a1bSJed Brown         ierr = PetscPrintf(PETSC_COMM_SELF," bus %D: load\n",i);CHKERRQ(ierr);
29*c4762a1bSJed Brown       }
30*c4762a1bSJed Brown     }
31*c4762a1bSJed Brown   }
32*c4762a1bSJed Brown   PetscFunctionReturn(0);
33*c4762a1bSJed Brown }
34*c4762a1bSJed Brown 
35*c4762a1bSJed Brown PetscErrorCode FormJacobian_Power_private(DM networkdm,Vec localX,Mat J,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx)
36*c4762a1bSJed Brown {
37*c4762a1bSJed Brown   PetscErrorCode    ierr;
38*c4762a1bSJed Brown   const PetscScalar *xarr;
39*c4762a1bSJed Brown   PetscInt          i,v,row[2],col[8],e,vfrom,vto;
40*c4762a1bSJed Brown   PetscInt          offsetfrom,offsetto,goffsetfrom,goffsetto,numComps;
41*c4762a1bSJed Brown   PetscScalar       values[8];
42*c4762a1bSJed Brown   PetscInt          j,key,offset,goffset;
43*c4762a1bSJed Brown   PetscScalar       Vm;
44*c4762a1bSJed Brown   UserCtx_Power     *user_power=(UserCtx_Power*)appctx;
45*c4762a1bSJed Brown   PetscScalar       Sbase=user_power->Sbase;
46*c4762a1bSJed Brown   VERTEX_Power      bus;
47*c4762a1bSJed Brown   PetscBool         ghostvtex;
48*c4762a1bSJed Brown   void*             component;
49*c4762a1bSJed Brown 
50*c4762a1bSJed Brown   PetscFunctionBegin;
51*c4762a1bSJed Brown   ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr);
52*c4762a1bSJed Brown 
53*c4762a1bSJed Brown   for (v=0; v<nv; v++) {
54*c4762a1bSJed Brown     ierr = DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex);CHKERRQ(ierr);
55*c4762a1bSJed Brown 
56*c4762a1bSJed Brown     ierr = DMNetworkGetNumComponents(networkdm,vtx[v],&numComps);CHKERRQ(ierr);
57*c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
58*c4762a1bSJed Brown       ierr = DMNetworkGetVariableOffset(networkdm,vtx[v],&offset);CHKERRQ(ierr);
59*c4762a1bSJed Brown       ierr = DMNetworkGetVariableGlobalOffset(networkdm,vtx[v],&goffset);CHKERRQ(ierr);
60*c4762a1bSJed Brown       ierr = DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component);CHKERRQ(ierr);
61*c4762a1bSJed Brown 
62*c4762a1bSJed Brown       if (key == user_power->compkey_bus) {
63*c4762a1bSJed Brown         PetscInt       nconnedges;
64*c4762a1bSJed Brown 	const PetscInt *connedges;
65*c4762a1bSJed Brown 
66*c4762a1bSJed Brown 	bus = (VERTEX_Power)(component);
67*c4762a1bSJed Brown         if (!ghostvtex) {
68*c4762a1bSJed Brown           /* Handle reference bus constrained dofs */
69*c4762a1bSJed Brown           if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
70*c4762a1bSJed Brown             row[0] = goffset; row[1] = goffset+1;
71*c4762a1bSJed Brown             col[0] = goffset; col[1] = goffset+1; col[2] = goffset; col[3] = goffset+1;
72*c4762a1bSJed Brown             values[0] = 1.0; values[1] = 0.0; values[2] = 0.0; values[3] = 1.0;
73*c4762a1bSJed Brown             ierr = MatSetValues(J,2,row,2,col,values,ADD_VALUES);CHKERRQ(ierr);
74*c4762a1bSJed Brown             break;
75*c4762a1bSJed Brown           }
76*c4762a1bSJed Brown 
77*c4762a1bSJed Brown           Vm = xarr[offset+1];
78*c4762a1bSJed Brown 
79*c4762a1bSJed Brown           /* Shunt injections */
80*c4762a1bSJed Brown           row[0] = goffset; row[1] = goffset+1;
81*c4762a1bSJed Brown           col[0] = goffset; col[1] = goffset+1;
82*c4762a1bSJed Brown           values[0] = values[1] = values[2] = values[3] = 0.0;
83*c4762a1bSJed Brown           if (bus->ide != PV_BUS) {
84*c4762a1bSJed Brown             values[1] = 2.0*Vm*bus->gl/Sbase;
85*c4762a1bSJed Brown             values[3] = -2.0*Vm*bus->bl/Sbase;
86*c4762a1bSJed Brown           }
87*c4762a1bSJed Brown           ierr = MatSetValues(J,2,row,2,col,values,ADD_VALUES);CHKERRQ(ierr);
88*c4762a1bSJed Brown         }
89*c4762a1bSJed Brown 
90*c4762a1bSJed Brown 	ierr = DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown 	for (i=0; i < nconnedges; i++) {
93*c4762a1bSJed Brown 	  EDGE_Power     branch;
94*c4762a1bSJed Brown 	  VERTEX_Power   busf,bust;
95*c4762a1bSJed Brown 	  PetscInt       keyf,keyt;
96*c4762a1bSJed Brown           PetscScalar    Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt;
97*c4762a1bSJed Brown           const PetscInt *cone;
98*c4762a1bSJed Brown           PetscScalar    Vmf,Vmt,thetaf,thetat,thetaft,thetatf;
99*c4762a1bSJed Brown 
100*c4762a1bSJed Brown 	  e = connedges[i];
101*c4762a1bSJed Brown 	  ierr = DMNetworkGetComponent(networkdm,e,0,&key,(void**)&branch);CHKERRQ(ierr);
102*c4762a1bSJed Brown 	  if (!branch->status) continue;
103*c4762a1bSJed Brown 
104*c4762a1bSJed Brown 	  Gff = branch->yff[0];
105*c4762a1bSJed Brown 	  Bff = branch->yff[1];
106*c4762a1bSJed Brown 	  Gft = branch->yft[0];
107*c4762a1bSJed Brown 	  Bft = branch->yft[1];
108*c4762a1bSJed Brown 	  Gtf = branch->ytf[0];
109*c4762a1bSJed Brown 	  Btf = branch->ytf[1];
110*c4762a1bSJed Brown 	  Gtt = branch->ytt[0];
111*c4762a1bSJed Brown 	  Btt = branch->ytt[1];
112*c4762a1bSJed Brown 
113*c4762a1bSJed Brown 	  ierr = DMNetworkGetConnectedVertices(networkdm,edges[e],&cone);CHKERRQ(ierr);
114*c4762a1bSJed Brown 	  vfrom = cone[0];
115*c4762a1bSJed Brown 	  vto   = cone[1];
116*c4762a1bSJed Brown 
117*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);CHKERRQ(ierr);
118*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableOffset(networkdm,vto,&offsetto);CHKERRQ(ierr);
119*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableGlobalOffset(networkdm,vfrom,&goffsetfrom);CHKERRQ(ierr);
120*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableGlobalOffset(networkdm,vto,&goffsetto);CHKERRQ(ierr);
121*c4762a1bSJed Brown 
122*c4762a1bSJed Brown 	  if (goffsetto < 0) goffsetto = -goffsetto - 1;
123*c4762a1bSJed Brown 
124*c4762a1bSJed Brown 	  thetaf  = xarr[offsetfrom];
125*c4762a1bSJed Brown 	  Vmf     = xarr[offsetfrom+1];
126*c4762a1bSJed Brown 	  thetat  = xarr[offsetto];
127*c4762a1bSJed Brown 	  Vmt     = xarr[offsetto+1];
128*c4762a1bSJed Brown 	  thetaft = thetaf - thetat;
129*c4762a1bSJed Brown 	  thetatf = thetat - thetaf;
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown 	  ierr = DMNetworkGetComponent(networkdm,vfrom,0,&keyf,(void**)&busf);CHKERRQ(ierr);
132*c4762a1bSJed Brown 	  ierr = DMNetworkGetComponent(networkdm,vto,0,&keyt,(void**)&bust);CHKERRQ(ierr);
133*c4762a1bSJed Brown 
134*c4762a1bSJed Brown 	  if (vfrom == vtx[v]) {
135*c4762a1bSJed Brown 	    if (busf->ide != REF_BUS) {
136*c4762a1bSJed Brown 	      /*    farr[offsetfrom]   += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft));  */
137*c4762a1bSJed Brown 	      row[0]  = goffsetfrom;
138*c4762a1bSJed Brown 	      col[0]  = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1;
139*c4762a1bSJed Brown 	      values[0] =  Vmf*Vmt*(Gft*-PetscSinScalar(thetaft) + Bft*PetscCosScalar(thetaft)); /* df_dthetaf */
140*c4762a1bSJed Brown 	      values[1] =  2.0*Gff*Vmf + Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmf */
141*c4762a1bSJed Brown 	      values[2] =  Vmf*Vmt*(Gft*PetscSinScalar(thetaft) + Bft*-PetscCosScalar(thetaft)); /* df_dthetat */
142*c4762a1bSJed Brown 	      values[3] =  Vmf*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmt */
143*c4762a1bSJed Brown 
144*c4762a1bSJed Brown 	      ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr);
145*c4762a1bSJed Brown 	    }
146*c4762a1bSJed Brown 	    if (busf->ide != PV_BUS && busf->ide != REF_BUS) {
147*c4762a1bSJed Brown 	      row[0] = goffsetfrom+1;
148*c4762a1bSJed Brown 	      col[0]  = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1;
149*c4762a1bSJed Brown 	      /*    farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */
150*c4762a1bSJed Brown 	      values[0] =  Vmf*Vmt*(Bft*PetscSinScalar(thetaft) + Gft*PetscCosScalar(thetaft));
151*c4762a1bSJed Brown 	      values[1] =  -2.0*Bff*Vmf + Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
152*c4762a1bSJed Brown 	      values[2] =  Vmf*Vmt*(-Bft*PetscSinScalar(thetaft) + Gft*-PetscCosScalar(thetaft));
153*c4762a1bSJed Brown 	      values[3] =  Vmf*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
154*c4762a1bSJed Brown 
155*c4762a1bSJed Brown 	      ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr);
156*c4762a1bSJed Brown 	    }
157*c4762a1bSJed Brown 	  } else {
158*c4762a1bSJed Brown 	    if (bust->ide != REF_BUS) {
159*c4762a1bSJed Brown 	      row[0] = goffsetto;
160*c4762a1bSJed Brown 	      col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1;
161*c4762a1bSJed Brown 	      /*    farr[offsetto]   += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */
162*c4762a1bSJed Brown 	      values[0] =  Vmt*Vmf*(Gtf*-PetscSinScalar(thetatf) + Btf*PetscCosScalar(thetaft)); /* df_dthetat */
163*c4762a1bSJed Brown 	      values[1] =  2.0*Gtt*Vmt + Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmt */
164*c4762a1bSJed Brown 	      values[2] =  Vmt*Vmf*(Gtf*PetscSinScalar(thetatf) + Btf*-PetscCosScalar(thetatf)); /* df_dthetaf */
165*c4762a1bSJed Brown 	      values[3] =  Vmt*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmf */
166*c4762a1bSJed Brown 
167*c4762a1bSJed Brown 	      ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr);
168*c4762a1bSJed Brown 	    }
169*c4762a1bSJed Brown 	    if (bust->ide != PV_BUS && bust->ide != REF_BUS) {
170*c4762a1bSJed Brown 	      row[0] = goffsetto+1;
171*c4762a1bSJed Brown 	      col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1;
172*c4762a1bSJed Brown 	      /*    farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */
173*c4762a1bSJed Brown 	      values[0] =  Vmt*Vmf*(Btf*PetscSinScalar(thetatf) + Gtf*PetscCosScalar(thetatf));
174*c4762a1bSJed Brown 	      values[1] =  -2.0*Btt*Vmt + Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
175*c4762a1bSJed Brown 	      values[2] =  Vmt*Vmf*(-Btf*PetscSinScalar(thetatf) + Gtf*-PetscCosScalar(thetatf));
176*c4762a1bSJed Brown 	      values[3] =  Vmt*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown 	      ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr);
179*c4762a1bSJed Brown 	    }
180*c4762a1bSJed Brown 	  }
181*c4762a1bSJed Brown 	}
182*c4762a1bSJed Brown 	if (!ghostvtex && bus->ide == PV_BUS) {
183*c4762a1bSJed Brown 	  row[0] = goffset+1; col[0] = goffset+1; values[0]  = 1.0;
184*c4762a1bSJed Brown           if (user_power->jac_error) values[0] = 50.0;
185*c4762a1bSJed Brown           ierr = MatSetValues(J,1,row,1,col,values,ADD_VALUES);CHKERRQ(ierr);
186*c4762a1bSJed Brown 	}
187*c4762a1bSJed Brown       }
188*c4762a1bSJed Brown     }
189*c4762a1bSJed Brown   }
190*c4762a1bSJed Brown 
191*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr);
192*c4762a1bSJed Brown   PetscFunctionReturn(0);
193*c4762a1bSJed Brown }
194*c4762a1bSJed Brown 
195*c4762a1bSJed Brown PetscErrorCode FormJacobian_Power(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx)
196*c4762a1bSJed Brown {
197*c4762a1bSJed Brown   PetscErrorCode ierr;
198*c4762a1bSJed Brown   DM             networkdm;
199*c4762a1bSJed Brown   Vec            localX;
200*c4762a1bSJed Brown   PetscInt       nv,ne;
201*c4762a1bSJed Brown   const PetscInt *vtx,*edges;
202*c4762a1bSJed Brown 
203*c4762a1bSJed Brown   PetscFunctionBegin;
204*c4762a1bSJed Brown   ierr = MatZeroEntries(J);CHKERRQ(ierr);
205*c4762a1bSJed Brown 
206*c4762a1bSJed Brown   ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr);
207*c4762a1bSJed Brown   ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr);
208*c4762a1bSJed Brown 
209*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
210*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
211*c4762a1bSJed Brown 
212*c4762a1bSJed Brown   ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
213*c4762a1bSJed Brown   ierr = FormJacobian_Power_private(networkdm,localX,J,nv,ne,vtx,edges,appctx);CHKERRQ(ierr);
214*c4762a1bSJed Brown 
215*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr);
216*c4762a1bSJed Brown 
217*c4762a1bSJed Brown   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
218*c4762a1bSJed Brown   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
219*c4762a1bSJed Brown   PetscFunctionReturn(0);
220*c4762a1bSJed Brown }
221*c4762a1bSJed Brown 
222*c4762a1bSJed Brown PetscErrorCode FormFunction_Power(DM networkdm,Vec localX, Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx)
223*c4762a1bSJed Brown {
224*c4762a1bSJed Brown   PetscErrorCode    ierr;
225*c4762a1bSJed Brown   UserCtx_Power     *User=(UserCtx_Power*)appctx;
226*c4762a1bSJed Brown   PetscInt          e,v,vfrom,vto;
227*c4762a1bSJed Brown   const PetscScalar *xarr;
228*c4762a1bSJed Brown   PetscScalar       *farr;
229*c4762a1bSJed Brown   PetscInt          offsetfrom,offsetto,offset,i,j,key,numComps;
230*c4762a1bSJed Brown   PetscScalar       Vm;
231*c4762a1bSJed Brown   PetscScalar       Sbase=User->Sbase;
232*c4762a1bSJed Brown   VERTEX_Power      bus=NULL;
233*c4762a1bSJed Brown   GEN               gen;
234*c4762a1bSJed Brown   LOAD              load;
235*c4762a1bSJed Brown   PetscBool         ghostvtex;
236*c4762a1bSJed Brown   void*             component;
237*c4762a1bSJed Brown 
238*c4762a1bSJed Brown   PetscFunctionBegin;
239*c4762a1bSJed Brown   ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr);
240*c4762a1bSJed Brown   ierr = VecGetArray(localF,&farr);CHKERRQ(ierr);
241*c4762a1bSJed Brown 
242*c4762a1bSJed Brown   for (v=0; v<nv; v++) {
243*c4762a1bSJed Brown     ierr = DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex);CHKERRQ(ierr);
244*c4762a1bSJed Brown     ierr = DMNetworkGetNumComponents(networkdm,vtx[v],&numComps);CHKERRQ(ierr);
245*c4762a1bSJed Brown     ierr = DMNetworkGetVariableOffset(networkdm,vtx[v],&offset);CHKERRQ(ierr);
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
248*c4762a1bSJed Brown       ierr = DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component);CHKERRQ(ierr);
249*c4762a1bSJed Brown       if (key == User->compkey_bus) {
250*c4762a1bSJed Brown         PetscInt       nconnedges;
251*c4762a1bSJed Brown 	const PetscInt *connedges;
252*c4762a1bSJed Brown 
253*c4762a1bSJed Brown 	bus = (VERTEX_Power)(component);
254*c4762a1bSJed Brown 	/* Handle reference bus constrained dofs */
255*c4762a1bSJed Brown 	if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
256*c4762a1bSJed Brown 	  farr[offset] = xarr[offset] - bus->va*PETSC_PI/180.0;
257*c4762a1bSJed Brown 	  farr[offset+1] = xarr[offset+1] - bus->vm;
258*c4762a1bSJed Brown 	  break;
259*c4762a1bSJed Brown 	}
260*c4762a1bSJed Brown 
261*c4762a1bSJed Brown 	if (!ghostvtex) {
262*c4762a1bSJed Brown 	  Vm = xarr[offset+1];
263*c4762a1bSJed Brown 
264*c4762a1bSJed Brown 	  /* Shunt injections */
265*c4762a1bSJed Brown 	  farr[offset] += Vm*Vm*bus->gl/Sbase;
266*c4762a1bSJed Brown 	  if(bus->ide != PV_BUS) farr[offset+1] += -Vm*Vm*bus->bl/Sbase;
267*c4762a1bSJed Brown 	}
268*c4762a1bSJed Brown 
269*c4762a1bSJed Brown 	ierr = DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges);CHKERRQ(ierr);
270*c4762a1bSJed Brown 	for (i=0; i < nconnedges; i++) {
271*c4762a1bSJed Brown 	  EDGE_Power       branch;
272*c4762a1bSJed Brown 	  PetscInt       keye;
273*c4762a1bSJed Brown           PetscScalar    Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt;
274*c4762a1bSJed Brown           const PetscInt *cone;
275*c4762a1bSJed Brown           PetscScalar    Vmf,Vmt,thetaf,thetat,thetaft,thetatf;
276*c4762a1bSJed Brown 
277*c4762a1bSJed Brown 	  e = connedges[i];
278*c4762a1bSJed Brown 	  ierr = DMNetworkGetComponent(networkdm,e,0,&keye,(void**)&branch);CHKERRQ(ierr);
279*c4762a1bSJed Brown 	  if (!branch->status) continue;
280*c4762a1bSJed Brown 	  Gff = branch->yff[0];
281*c4762a1bSJed Brown 	  Bff = branch->yff[1];
282*c4762a1bSJed Brown 	  Gft = branch->yft[0];
283*c4762a1bSJed Brown 	  Bft = branch->yft[1];
284*c4762a1bSJed Brown 	  Gtf = branch->ytf[0];
285*c4762a1bSJed Brown 	  Btf = branch->ytf[1];
286*c4762a1bSJed Brown 	  Gtt = branch->ytt[0];
287*c4762a1bSJed Brown 	  Btt = branch->ytt[1];
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown 	  ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr);
290*c4762a1bSJed Brown 	  vfrom = cone[0];
291*c4762a1bSJed Brown 	  vto   = cone[1];
292*c4762a1bSJed Brown 
293*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);CHKERRQ(ierr);
294*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableOffset(networkdm,vto,&offsetto);CHKERRQ(ierr);
295*c4762a1bSJed Brown 
296*c4762a1bSJed Brown 	  thetaf = xarr[offsetfrom];
297*c4762a1bSJed Brown 	  Vmf     = xarr[offsetfrom+1];
298*c4762a1bSJed Brown 	  thetat = xarr[offsetto];
299*c4762a1bSJed Brown 	  Vmt     = xarr[offsetto+1];
300*c4762a1bSJed Brown 	  thetaft = thetaf - thetat;
301*c4762a1bSJed Brown 	  thetatf = thetat - thetaf;
302*c4762a1bSJed Brown 
303*c4762a1bSJed Brown 	  if (vfrom == vtx[v]) {
304*c4762a1bSJed Brown 	    farr[offsetfrom]   += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft));
305*c4762a1bSJed Brown 	    farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
306*c4762a1bSJed Brown 	  } else {
307*c4762a1bSJed Brown 	    farr[offsetto]   += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf));
308*c4762a1bSJed Brown 	    farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
309*c4762a1bSJed Brown 	  }
310*c4762a1bSJed Brown 	}
311*c4762a1bSJed Brown       } else if (key == User->compkey_gen) {
312*c4762a1bSJed Brown 	if (!ghostvtex) {
313*c4762a1bSJed Brown 	  gen = (GEN)(component);
314*c4762a1bSJed Brown 	  if (!gen->status) continue;
315*c4762a1bSJed Brown 	  farr[offset] += -gen->pg/Sbase;
316*c4762a1bSJed Brown 	  farr[offset+1] += -gen->qg/Sbase;
317*c4762a1bSJed Brown 	}
318*c4762a1bSJed Brown       } else if (key == User->compkey_load) {
319*c4762a1bSJed Brown 	if (!ghostvtex) {
320*c4762a1bSJed Brown 	  load = (LOAD)(component);
321*c4762a1bSJed Brown 	  farr[offset] += load->pl/Sbase;
322*c4762a1bSJed Brown 	  farr[offset+1] += load->ql/Sbase;
323*c4762a1bSJed Brown 	}
324*c4762a1bSJed Brown       }
325*c4762a1bSJed Brown     }
326*c4762a1bSJed Brown     if (bus && bus->ide == PV_BUS) {
327*c4762a1bSJed Brown       farr[offset+1] = xarr[offset+1] - bus->vm;
328*c4762a1bSJed Brown     }
329*c4762a1bSJed Brown   }
330*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr);
331*c4762a1bSJed Brown   ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr);
332*c4762a1bSJed Brown   PetscFunctionReturn(0);
333*c4762a1bSJed Brown }
334*c4762a1bSJed Brown 
335*c4762a1bSJed Brown PetscErrorCode SetInitialGuess_Power(DM networkdm,Vec localX,PetscInt nv,PetscInt ne,const PetscInt *vtx,const PetscInt *edges,void* appctx)
336*c4762a1bSJed Brown {
337*c4762a1bSJed Brown   PetscErrorCode ierr;
338*c4762a1bSJed Brown   VERTEX_Power   bus;
339*c4762a1bSJed Brown   PetscInt       i;
340*c4762a1bSJed Brown   GEN            gen;
341*c4762a1bSJed Brown   PetscBool      ghostvtex;
342*c4762a1bSJed Brown   PetscScalar    *xarr;
343*c4762a1bSJed Brown   PetscInt       key,numComps,j,offset;
344*c4762a1bSJed Brown   void*          component;
345*c4762a1bSJed Brown   PetscMPIInt    rank;
346*c4762a1bSJed Brown   MPI_Comm       comm;
347*c4762a1bSJed Brown   UserCtx_Power  *User=(UserCtx_Power*)appctx;
348*c4762a1bSJed Brown 
349*c4762a1bSJed Brown   PetscFunctionBegin;
350*c4762a1bSJed Brown   ierr = PetscObjectGetComm((PetscObject)networkdm,&comm);CHKERRQ(ierr);
351*c4762a1bSJed Brown   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
352*c4762a1bSJed Brown   ierr = VecGetArray(localX,&xarr);CHKERRQ(ierr);
353*c4762a1bSJed Brown   for (i = 0; i < nv; i++) {
354*c4762a1bSJed Brown     ierr = DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);CHKERRQ(ierr);
355*c4762a1bSJed Brown     if (ghostvtex) continue;
356*c4762a1bSJed Brown 
357*c4762a1bSJed Brown     ierr = DMNetworkGetVariableOffset(networkdm,vtx[i],&offset);CHKERRQ(ierr);
358*c4762a1bSJed Brown     ierr = DMNetworkGetNumComponents(networkdm,vtx[i],&numComps);CHKERRQ(ierr);
359*c4762a1bSJed Brown     for (j=0; j < numComps; j++) {
360*c4762a1bSJed Brown       ierr = DMNetworkGetComponent(networkdm,vtx[i],j,&key,&component);CHKERRQ(ierr);
361*c4762a1bSJed Brown       if (key == User->compkey_bus) {
362*c4762a1bSJed Brown 	bus = (VERTEX_Power)(component);
363*c4762a1bSJed Brown 	xarr[offset] = bus->va*PETSC_PI/180.0;
364*c4762a1bSJed Brown 	xarr[offset+1] = bus->vm;
365*c4762a1bSJed Brown       } else if(key == User->compkey_gen) {
366*c4762a1bSJed Brown 	gen = (GEN)(component);
367*c4762a1bSJed Brown 	if (!gen->status) continue;
368*c4762a1bSJed Brown 	xarr[offset+1] = gen->vs;
369*c4762a1bSJed Brown 	break;
370*c4762a1bSJed Brown       }
371*c4762a1bSJed Brown     }
372*c4762a1bSJed Brown   }
373*c4762a1bSJed Brown   ierr = VecRestoreArray(localX,&xarr);CHKERRQ(ierr);
374*c4762a1bSJed Brown   PetscFunctionReturn(0);
375*c4762a1bSJed Brown }
376