xref: /petsc/src/snes/tutorials/network/power/power2.c (revision c4762a1b19cd2af06abeed90e8f9d34fb975dd94)
1*c4762a1bSJed Brown static char help[] = "This example demonstrates the use of DMNetwork interface with subnetworks for solving a nonlinear electric power grid problem.\n\
2*c4762a1bSJed Brown                       The available solver options are in the poweroptions file and the data files are in the datafiles directory.\n\
3*c4762a1bSJed Brown                       The data file format used is from the MatPower package (http://www.pserc.cornell.edu//matpower/).\n\
4*c4762a1bSJed Brown                       This example shows the use of subnetwork feature in DMNetwork. It creates duplicates of the same network which are treated as subnetworks.\n\
5*c4762a1bSJed Brown                       Run this program: mpiexec -n <n> ./pf2\n\
6*c4762a1bSJed Brown                       mpiexec -n <n> ./pf2 \n";
7*c4762a1bSJed Brown 
8*c4762a1bSJed Brown /* T
9*c4762a1bSJed Brown    Concepts: DMNetwork
10*c4762a1bSJed Brown    Concepts: PETSc SNES solver
11*c4762a1bSJed Brown */
12*c4762a1bSJed Brown 
13*c4762a1bSJed Brown #include "power.h"
14*c4762a1bSJed Brown #include <petscdmnetwork.h>
15*c4762a1bSJed Brown 
16*c4762a1bSJed Brown PetscErrorCode FormFunction_Subnet(DM networkdm,Vec localX, Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx)
17*c4762a1bSJed Brown {
18*c4762a1bSJed Brown   PetscErrorCode ierr;
19*c4762a1bSJed Brown   UserCtx_Power  *User=(UserCtx_Power*)appctx;
20*c4762a1bSJed Brown   PetscInt       e,v,vfrom,vto;
21*c4762a1bSJed Brown   const PetscScalar *xarr;
22*c4762a1bSJed Brown   PetscScalar    *farr;
23*c4762a1bSJed Brown   PetscInt       offsetfrom,offsetto,offset;
24*c4762a1bSJed Brown 
25*c4762a1bSJed Brown   PetscFunctionBegin;
26*c4762a1bSJed Brown   ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr);
27*c4762a1bSJed Brown   ierr = VecGetArray(localF,&farr);CHKERRQ(ierr);
28*c4762a1bSJed Brown 
29*c4762a1bSJed Brown   for (v=0; v<nv; v++) {
30*c4762a1bSJed Brown     PetscInt    i,j,key;
31*c4762a1bSJed Brown     PetscScalar Vm;
32*c4762a1bSJed Brown     PetscScalar Sbase=User->Sbase;
33*c4762a1bSJed Brown     VERTEX_Power  bus=NULL;
34*c4762a1bSJed Brown     GEN         gen;
35*c4762a1bSJed Brown     LOAD        load;
36*c4762a1bSJed Brown     PetscBool   ghostvtex;
37*c4762a1bSJed Brown     PetscInt    numComps;
38*c4762a1bSJed Brown     void*       component;
39*c4762a1bSJed Brown 
40*c4762a1bSJed Brown     ierr = DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex);CHKERRQ(ierr);
41*c4762a1bSJed Brown     ierr = DMNetworkGetNumComponents(networkdm,vtx[v],&numComps);CHKERRQ(ierr);
42*c4762a1bSJed Brown     ierr = DMNetworkGetVariableOffset(networkdm,vtx[v],&offset);CHKERRQ(ierr);
43*c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
44*c4762a1bSJed Brown       ierr = DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component);CHKERRQ(ierr);
45*c4762a1bSJed Brown       if (key == 1) {
46*c4762a1bSJed Brown         PetscInt       nconnedges;
47*c4762a1bSJed Brown 	const PetscInt *connedges;
48*c4762a1bSJed Brown 
49*c4762a1bSJed Brown 	bus = (VERTEX_Power)(component);
50*c4762a1bSJed Brown 	/* Handle reference bus constrained dofs */
51*c4762a1bSJed Brown 	if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
52*c4762a1bSJed Brown 	  farr[offset] = xarr[offset] - bus->va*PETSC_PI/180.0;
53*c4762a1bSJed Brown 	  farr[offset+1] = xarr[offset+1] - bus->vm;
54*c4762a1bSJed Brown 	  break;
55*c4762a1bSJed Brown 	}
56*c4762a1bSJed Brown 
57*c4762a1bSJed Brown 	if (!ghostvtex) {
58*c4762a1bSJed Brown 	  Vm = xarr[offset+1];
59*c4762a1bSJed Brown 
60*c4762a1bSJed Brown 	  /* Shunt injections */
61*c4762a1bSJed Brown 	  farr[offset] += Vm*Vm*bus->gl/Sbase;
62*c4762a1bSJed Brown 	  if(bus->ide != PV_BUS) farr[offset+1] += -Vm*Vm*bus->bl/Sbase;
63*c4762a1bSJed Brown 	}
64*c4762a1bSJed Brown 
65*c4762a1bSJed Brown 	ierr = DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges);CHKERRQ(ierr);
66*c4762a1bSJed Brown 	for (i=0; i < nconnedges; i++) {
67*c4762a1bSJed Brown 	  EDGE_Power       branch;
68*c4762a1bSJed Brown 	  PetscInt       keye;
69*c4762a1bSJed Brown           PetscScalar    Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt;
70*c4762a1bSJed Brown           const PetscInt *cone;
71*c4762a1bSJed Brown           PetscScalar    Vmf,Vmt,thetaf,thetat,thetaft,thetatf;
72*c4762a1bSJed Brown 
73*c4762a1bSJed Brown 	  e = connedges[i];
74*c4762a1bSJed Brown 	  ierr = DMNetworkGetComponent(networkdm,e,0,&keye,(void**)&branch);CHKERRQ(ierr);
75*c4762a1bSJed Brown 	  if (!branch->status) continue;
76*c4762a1bSJed Brown 	  Gff = branch->yff[0];
77*c4762a1bSJed Brown 	  Bff = branch->yff[1];
78*c4762a1bSJed Brown 	  Gft = branch->yft[0];
79*c4762a1bSJed Brown 	  Bft = branch->yft[1];
80*c4762a1bSJed Brown 	  Gtf = branch->ytf[0];
81*c4762a1bSJed Brown 	  Btf = branch->ytf[1];
82*c4762a1bSJed Brown 	  Gtt = branch->ytt[0];
83*c4762a1bSJed Brown 	  Btt = branch->ytt[1];
84*c4762a1bSJed Brown 
85*c4762a1bSJed Brown 	  ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr);
86*c4762a1bSJed Brown 	  vfrom = cone[0];
87*c4762a1bSJed Brown 	  vto   = cone[1];
88*c4762a1bSJed Brown 
89*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);CHKERRQ(ierr);
90*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableOffset(networkdm,vto,&offsetto);CHKERRQ(ierr);
91*c4762a1bSJed Brown 
92*c4762a1bSJed Brown 	  thetaf = xarr[offsetfrom];
93*c4762a1bSJed Brown 	  Vmf     = xarr[offsetfrom+1];
94*c4762a1bSJed Brown 	  thetat = xarr[offsetto];
95*c4762a1bSJed Brown 	  Vmt     = xarr[offsetto+1];
96*c4762a1bSJed Brown 	  thetaft = thetaf - thetat;
97*c4762a1bSJed Brown 	  thetatf = thetat - thetaf;
98*c4762a1bSJed Brown 
99*c4762a1bSJed Brown 	  if (vfrom == vtx[v]) {
100*c4762a1bSJed Brown 	    farr[offsetfrom]   += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft));
101*c4762a1bSJed Brown 	    farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
102*c4762a1bSJed Brown 	  } else {
103*c4762a1bSJed Brown 	    farr[offsetto]   += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf));
104*c4762a1bSJed Brown 	    farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
105*c4762a1bSJed Brown 	  }
106*c4762a1bSJed Brown 	}
107*c4762a1bSJed Brown       } else if (key == 2) {
108*c4762a1bSJed Brown 	if (!ghostvtex) {
109*c4762a1bSJed Brown 	  gen = (GEN)(component);
110*c4762a1bSJed Brown 	  if (!gen->status) continue;
111*c4762a1bSJed Brown 	  farr[offset] += -gen->pg/Sbase;
112*c4762a1bSJed Brown 	  farr[offset+1] += -gen->qg/Sbase;
113*c4762a1bSJed Brown 	}
114*c4762a1bSJed Brown       } else if (key == 3) {
115*c4762a1bSJed Brown 	if (!ghostvtex) {
116*c4762a1bSJed Brown 	  load = (LOAD)(component);
117*c4762a1bSJed Brown 	  farr[offset] += load->pl/Sbase;
118*c4762a1bSJed Brown 	  farr[offset+1] += load->ql/Sbase;
119*c4762a1bSJed Brown 	}
120*c4762a1bSJed Brown       }
121*c4762a1bSJed Brown     }
122*c4762a1bSJed Brown     if (bus && bus->ide == PV_BUS) {
123*c4762a1bSJed Brown       farr[offset+1] = xarr[offset+1] - bus->vm;
124*c4762a1bSJed Brown     }
125*c4762a1bSJed Brown   }
126*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr);
127*c4762a1bSJed Brown   ierr = VecRestoreArray(localF,&farr);CHKERRQ(ierr);
128*c4762a1bSJed Brown   PetscFunctionReturn(0);
129*c4762a1bSJed Brown }
130*c4762a1bSJed Brown 
131*c4762a1bSJed Brown 
132*c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X, Vec F,void *appctx)
133*c4762a1bSJed Brown {
134*c4762a1bSJed Brown   PetscErrorCode ierr;
135*c4762a1bSJed Brown   DM             networkdm;
136*c4762a1bSJed Brown   Vec           localX,localF;
137*c4762a1bSJed Brown   PetscInt      nv,ne;
138*c4762a1bSJed Brown   const PetscInt *vtx,*edges;
139*c4762a1bSJed Brown 
140*c4762a1bSJed Brown   PetscFunctionBegin;
141*c4762a1bSJed Brown   ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr);
142*c4762a1bSJed Brown   ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr);
143*c4762a1bSJed Brown   ierr = DMGetLocalVector(networkdm,&localF);CHKERRQ(ierr);
144*c4762a1bSJed Brown   ierr = VecSet(F,0.0);CHKERRQ(ierr);
145*c4762a1bSJed Brown 
146*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
147*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
148*c4762a1bSJed Brown 
149*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(networkdm,F,INSERT_VALUES,localF);CHKERRQ(ierr);
150*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(networkdm,F,INSERT_VALUES,localF);CHKERRQ(ierr);
151*c4762a1bSJed Brown 
152*c4762a1bSJed Brown   /* Form Function for first subnetwork */
153*c4762a1bSJed Brown   ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
154*c4762a1bSJed Brown   ierr = FormFunction_Subnet(networkdm,localX,localF,nv,ne,vtx,edges,appctx);CHKERRQ(ierr);
155*c4762a1bSJed Brown 
156*c4762a1bSJed Brown   /* Form Function for second subnetwork */
157*c4762a1bSJed Brown   ierr = DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
158*c4762a1bSJed Brown   ierr = FormFunction_Subnet(networkdm,localX,localF,nv,ne,vtx,edges,appctx);CHKERRQ(ierr);
159*c4762a1bSJed Brown 
160*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr);
161*c4762a1bSJed Brown 
162*c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr);
163*c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F);CHKERRQ(ierr);
164*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(networkdm,&localF);CHKERRQ(ierr);
165*c4762a1bSJed Brown   PetscFunctionReturn(0);
166*c4762a1bSJed Brown }
167*c4762a1bSJed Brown 
168*c4762a1bSJed Brown PetscErrorCode FormJacobian_Subnet(DM networkdm,Vec localX, Mat J, Mat Jpre, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
169*c4762a1bSJed Brown {
170*c4762a1bSJed Brown   PetscErrorCode ierr;
171*c4762a1bSJed Brown   UserCtx_Power  *User=(UserCtx_Power*)appctx;
172*c4762a1bSJed Brown   PetscInt       e,v,vfrom,vto;
173*c4762a1bSJed Brown   const PetscScalar *xarr;
174*c4762a1bSJed Brown   PetscInt       offsetfrom,offsetto,goffsetfrom,goffsetto;
175*c4762a1bSJed Brown   PetscInt       row[2],col[8];
176*c4762a1bSJed Brown   PetscScalar    values[8];
177*c4762a1bSJed Brown 
178*c4762a1bSJed Brown   PetscFunctionBegin;
179*c4762a1bSJed Brown   ierr = VecGetArrayRead(localX,&xarr);CHKERRQ(ierr);
180*c4762a1bSJed Brown 
181*c4762a1bSJed Brown   for (v=0; v<nv; v++) {
182*c4762a1bSJed Brown     PetscInt    i,j,key;
183*c4762a1bSJed Brown     PetscInt    offset,goffset;
184*c4762a1bSJed Brown     PetscScalar Vm;
185*c4762a1bSJed Brown     PetscScalar Sbase=User->Sbase;
186*c4762a1bSJed Brown     VERTEX_Power  bus;
187*c4762a1bSJed Brown     PetscBool   ghostvtex;
188*c4762a1bSJed Brown     PetscInt    numComps;
189*c4762a1bSJed Brown     void*       component;
190*c4762a1bSJed Brown 
191*c4762a1bSJed Brown     ierr = DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex);CHKERRQ(ierr);
192*c4762a1bSJed Brown     ierr = DMNetworkGetNumComponents(networkdm,vtx[v],&numComps);CHKERRQ(ierr);
193*c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
194*c4762a1bSJed Brown       ierr = DMNetworkGetVariableOffset(networkdm,vtx[v],&offset);CHKERRQ(ierr);
195*c4762a1bSJed Brown       ierr = DMNetworkGetVariableGlobalOffset(networkdm,vtx[v],&goffset);CHKERRQ(ierr);
196*c4762a1bSJed Brown       ierr = DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component);CHKERRQ(ierr);
197*c4762a1bSJed Brown       if (key == 1) {
198*c4762a1bSJed Brown         PetscInt       nconnedges;
199*c4762a1bSJed Brown 	const PetscInt *connedges;
200*c4762a1bSJed Brown 
201*c4762a1bSJed Brown 	bus = (VERTEX_Power)(component);
202*c4762a1bSJed Brown 	if (!ghostvtex) {
203*c4762a1bSJed Brown 	  /* Handle reference bus constrained dofs */
204*c4762a1bSJed Brown 	  if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
205*c4762a1bSJed Brown 	    row[0] = goffset; row[1] = goffset+1;
206*c4762a1bSJed Brown 	    col[0] = goffset; col[1] = goffset+1; col[2] = goffset; col[3] = goffset+1;
207*c4762a1bSJed Brown 	    values[0] = 1.0; values[1] = 0.0; values[2] = 0.0; values[3] = 1.0;
208*c4762a1bSJed Brown 	    ierr = MatSetValues(J,2,row,2,col,values,ADD_VALUES);CHKERRQ(ierr);
209*c4762a1bSJed Brown 	    break;
210*c4762a1bSJed Brown 	  }
211*c4762a1bSJed Brown 
212*c4762a1bSJed Brown 	  Vm = xarr[offset+1];
213*c4762a1bSJed Brown 
214*c4762a1bSJed Brown 	  /* Shunt injections */
215*c4762a1bSJed Brown           row[0] = goffset; row[1] = goffset+1;
216*c4762a1bSJed Brown           col[0] = goffset; col[1] = goffset+1;
217*c4762a1bSJed Brown           values[0] = values[1] = values[2] = values[3] = 0.0;
218*c4762a1bSJed Brown           if (bus->ide != PV_BUS) {
219*c4762a1bSJed Brown             values[1] = 2.0*Vm*bus->gl/Sbase;
220*c4762a1bSJed Brown             values[3] = -2.0*Vm*bus->bl/Sbase;
221*c4762a1bSJed Brown           }
222*c4762a1bSJed Brown           ierr = MatSetValues(J,2,row,2,col,values,ADD_VALUES);CHKERRQ(ierr);
223*c4762a1bSJed Brown 	}
224*c4762a1bSJed Brown 
225*c4762a1bSJed Brown 	ierr = DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges);CHKERRQ(ierr);
226*c4762a1bSJed Brown 	for (i=0; i < nconnedges; i++) {
227*c4762a1bSJed Brown 	  EDGE_Power       branch;
228*c4762a1bSJed Brown 	  VERTEX_Power     busf,bust;
229*c4762a1bSJed Brown 	  PetscInt       keyf,keyt;
230*c4762a1bSJed Brown           PetscScalar    Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt;
231*c4762a1bSJed Brown           const PetscInt *cone;
232*c4762a1bSJed Brown           PetscScalar    Vmf,Vmt,thetaf,thetat,thetaft,thetatf;
233*c4762a1bSJed Brown 
234*c4762a1bSJed Brown 	  e = connedges[i];
235*c4762a1bSJed Brown 	  ierr = DMNetworkGetComponent(networkdm,e,0,&key,(void**)&branch);CHKERRQ(ierr);
236*c4762a1bSJed Brown 	  if (!branch->status) continue;
237*c4762a1bSJed Brown 
238*c4762a1bSJed Brown 	  Gff = branch->yff[0];
239*c4762a1bSJed Brown 	  Bff = branch->yff[1];
240*c4762a1bSJed Brown 	  Gft = branch->yft[0];
241*c4762a1bSJed Brown 	  Bft = branch->yft[1];
242*c4762a1bSJed Brown 	  Gtf = branch->ytf[0];
243*c4762a1bSJed Brown 	  Btf = branch->ytf[1];
244*c4762a1bSJed Brown 	  Gtt = branch->ytt[0];
245*c4762a1bSJed Brown 	  Btt = branch->ytt[1];
246*c4762a1bSJed Brown 
247*c4762a1bSJed Brown 	  ierr = DMNetworkGetConnectedVertices(networkdm,e,&cone);CHKERRQ(ierr);
248*c4762a1bSJed Brown 	  vfrom = cone[0];
249*c4762a1bSJed Brown 	  vto   = cone[1];
250*c4762a1bSJed Brown 
251*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableOffset(networkdm,vfrom,&offsetfrom);CHKERRQ(ierr);
252*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableOffset(networkdm,vto,&offsetto);CHKERRQ(ierr);
253*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableGlobalOffset(networkdm,vfrom,&goffsetfrom);CHKERRQ(ierr);
254*c4762a1bSJed Brown 	  ierr = DMNetworkGetVariableGlobalOffset(networkdm,vto,&goffsetto);CHKERRQ(ierr);
255*c4762a1bSJed Brown 
256*c4762a1bSJed Brown 	  if (goffsetto < 0) goffsetto = -goffsetto - 1;
257*c4762a1bSJed Brown 
258*c4762a1bSJed Brown 	  thetaf = xarr[offsetfrom];
259*c4762a1bSJed Brown 	  Vmf     = xarr[offsetfrom+1];
260*c4762a1bSJed Brown 	  thetat = xarr[offsetto];
261*c4762a1bSJed Brown 	  Vmt     = xarr[offsetto+1];
262*c4762a1bSJed Brown 	  thetaft = thetaf - thetat;
263*c4762a1bSJed Brown 	  thetatf = thetat - thetaf;
264*c4762a1bSJed Brown 
265*c4762a1bSJed Brown 	  ierr = DMNetworkGetComponent(networkdm,vfrom,0,&keyf,(void**)&busf);CHKERRQ(ierr);
266*c4762a1bSJed Brown 	  ierr = DMNetworkGetComponent(networkdm,vto,0,&keyt,(void**)&bust);CHKERRQ(ierr);
267*c4762a1bSJed Brown 
268*c4762a1bSJed Brown 	  if (vfrom == vtx[v]) {
269*c4762a1bSJed Brown 	    if (busf->ide != REF_BUS) {
270*c4762a1bSJed Brown 	      /*    farr[offsetfrom]   += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft));  */
271*c4762a1bSJed Brown 	      row[0]  = goffsetfrom;
272*c4762a1bSJed Brown 	      col[0]  = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1;
273*c4762a1bSJed Brown 	      values[0] =  Vmf*Vmt*(Gft*-PetscSinScalar(thetaft) + Bft*PetscCosScalar(thetaft)); /* df_dthetaf */
274*c4762a1bSJed Brown 	      values[1] =  2.0*Gff*Vmf + Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmf */
275*c4762a1bSJed Brown 	      values[2] =  Vmf*Vmt*(Gft*PetscSinScalar(thetaft) + Bft*-PetscCosScalar(thetaft)); /* df_dthetat */
276*c4762a1bSJed Brown 	      values[3] =  Vmf*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmt */
277*c4762a1bSJed Brown 
278*c4762a1bSJed Brown 	      ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr);
279*c4762a1bSJed Brown 	    }
280*c4762a1bSJed Brown 	    if (busf->ide != PV_BUS && busf->ide != REF_BUS) {
281*c4762a1bSJed Brown 	      row[0] = goffsetfrom+1;
282*c4762a1bSJed Brown 	      col[0]  = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1;
283*c4762a1bSJed Brown 	      /*    farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */
284*c4762a1bSJed Brown 	      values[0] =  Vmf*Vmt*(Bft*PetscSinScalar(thetaft) + Gft*PetscCosScalar(thetaft));
285*c4762a1bSJed Brown 	      values[1] =  -2.0*Bff*Vmf + Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
286*c4762a1bSJed Brown 	      values[2] =  Vmf*Vmt*(-Bft*PetscSinScalar(thetaft) + Gft*-PetscCosScalar(thetaft));
287*c4762a1bSJed Brown 	      values[3] =  Vmf*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
288*c4762a1bSJed Brown 
289*c4762a1bSJed Brown 	      ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr);
290*c4762a1bSJed Brown 	    }
291*c4762a1bSJed Brown 	  } else {
292*c4762a1bSJed Brown 	    if (bust->ide != REF_BUS) {
293*c4762a1bSJed Brown 	      row[0] = goffsetto;
294*c4762a1bSJed Brown 	      col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1;
295*c4762a1bSJed Brown 	      /*    farr[offsetto]   += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */
296*c4762a1bSJed Brown 	      values[0] =  Vmt*Vmf*(Gtf*-PetscSinScalar(thetatf) + Btf*PetscCosScalar(thetaft)); /* df_dthetat */
297*c4762a1bSJed Brown 	      values[1] =  2.0*Gtt*Vmt + Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmt */
298*c4762a1bSJed Brown 	      values[2] =  Vmt*Vmf*(Gtf*PetscSinScalar(thetatf) + Btf*-PetscCosScalar(thetatf)); /* df_dthetaf */
299*c4762a1bSJed Brown 	      values[3] =  Vmt*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmf */
300*c4762a1bSJed Brown 
301*c4762a1bSJed Brown 	      ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr);
302*c4762a1bSJed Brown 	    }
303*c4762a1bSJed Brown 	    if (bust->ide != PV_BUS && bust->ide != REF_BUS) {
304*c4762a1bSJed Brown 	      row[0] = goffsetto+1;
305*c4762a1bSJed Brown 	      col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1;
306*c4762a1bSJed Brown 	      /*    farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */
307*c4762a1bSJed Brown 	      values[0] =  Vmt*Vmf*(Btf*PetscSinScalar(thetatf) + Gtf*PetscCosScalar(thetatf));
308*c4762a1bSJed Brown 	      values[1] =  -2.0*Btt*Vmt + Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
309*c4762a1bSJed Brown 	      values[2] =  Vmt*Vmf*(-Btf*PetscSinScalar(thetatf) + Gtf*-PetscCosScalar(thetatf));
310*c4762a1bSJed Brown 	      values[3] =  Vmt*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
311*c4762a1bSJed Brown 
312*c4762a1bSJed Brown 	      ierr = MatSetValues(J,1,row,4,col,values,ADD_VALUES);CHKERRQ(ierr);
313*c4762a1bSJed Brown 	    }
314*c4762a1bSJed Brown 	  }
315*c4762a1bSJed Brown 	}
316*c4762a1bSJed Brown 	if (!ghostvtex && bus->ide == PV_BUS) {
317*c4762a1bSJed Brown 	  row[0] = goffset+1; col[0] = goffset+1;
318*c4762a1bSJed Brown 	  values[0]  = 1.0;
319*c4762a1bSJed Brown 	  ierr = MatSetValues(J,1,row,1,col,values,ADD_VALUES);CHKERRQ(ierr);
320*c4762a1bSJed Brown 	}
321*c4762a1bSJed Brown       }
322*c4762a1bSJed Brown     }
323*c4762a1bSJed Brown   }
324*c4762a1bSJed Brown   ierr = VecRestoreArrayRead(localX,&xarr);CHKERRQ(ierr);
325*c4762a1bSJed Brown   PetscFunctionReturn(0);
326*c4762a1bSJed Brown }
327*c4762a1bSJed Brown 
328*c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx)
329*c4762a1bSJed Brown {
330*c4762a1bSJed Brown   PetscErrorCode ierr;
331*c4762a1bSJed Brown   DM            networkdm;
332*c4762a1bSJed Brown   Vec           localX;
333*c4762a1bSJed Brown   PetscInt      ne,nv;
334*c4762a1bSJed Brown   const PetscInt *vtx,*edges;
335*c4762a1bSJed Brown 
336*c4762a1bSJed Brown   PetscFunctionBegin;
337*c4762a1bSJed Brown   ierr = MatZeroEntries(J);CHKERRQ(ierr);
338*c4762a1bSJed Brown 
339*c4762a1bSJed Brown   ierr = SNESGetDM(snes,&networkdm);CHKERRQ(ierr);
340*c4762a1bSJed Brown   ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr);
341*c4762a1bSJed Brown 
342*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
343*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
344*c4762a1bSJed Brown 
345*c4762a1bSJed Brown   /* Form Jacobian for first subnetwork */
346*c4762a1bSJed Brown   ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
347*c4762a1bSJed Brown   ierr = FormJacobian_Subnet(networkdm,localX,J,Jpre,nv,ne,vtx,edges,appctx);CHKERRQ(ierr);
348*c4762a1bSJed Brown 
349*c4762a1bSJed Brown   /* Form Jacobian for second subnetwork */
350*c4762a1bSJed Brown   ierr = DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
351*c4762a1bSJed Brown   ierr = FormJacobian_Subnet(networkdm,localX,J,Jpre,nv,ne,vtx,edges,appctx);CHKERRQ(ierr);
352*c4762a1bSJed Brown 
353*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr);
354*c4762a1bSJed Brown 
355*c4762a1bSJed Brown   ierr = MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
356*c4762a1bSJed Brown   ierr = MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
357*c4762a1bSJed Brown   PetscFunctionReturn(0);
358*c4762a1bSJed Brown }
359*c4762a1bSJed Brown 
360*c4762a1bSJed Brown PetscErrorCode SetInitialValues_Subnet(DM networkdm,Vec localX,PetscInt nv,PetscInt ne, const PetscInt *vtx, const PetscInt *edges,void* appctx)
361*c4762a1bSJed Brown {
362*c4762a1bSJed Brown   PetscErrorCode ierr;
363*c4762a1bSJed Brown   VERTEX_Power     bus;
364*c4762a1bSJed Brown   PetscInt       i;
365*c4762a1bSJed Brown   GEN            gen;
366*c4762a1bSJed Brown   PetscBool      ghostvtex;
367*c4762a1bSJed Brown   PetscScalar    *xarr;
368*c4762a1bSJed Brown   PetscInt       key,numComps,j,offset;
369*c4762a1bSJed Brown   void*          component;
370*c4762a1bSJed Brown 
371*c4762a1bSJed Brown   PetscFunctionBegin;
372*c4762a1bSJed Brown   ierr = VecGetArray(localX,&xarr);CHKERRQ(ierr);
373*c4762a1bSJed Brown   for (i = 0; i < nv; i++) {
374*c4762a1bSJed Brown     ierr = DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex);CHKERRQ(ierr);
375*c4762a1bSJed Brown     if (ghostvtex) continue;
376*c4762a1bSJed Brown 
377*c4762a1bSJed Brown     ierr = DMNetworkGetVariableOffset(networkdm,vtx[i],&offset);CHKERRQ(ierr);
378*c4762a1bSJed Brown     ierr = DMNetworkGetNumComponents(networkdm,vtx[i],&numComps);CHKERRQ(ierr);
379*c4762a1bSJed Brown     for (j=0; j < numComps; j++) {
380*c4762a1bSJed Brown       ierr = DMNetworkGetComponent(networkdm,vtx[i],j,&key,&component);CHKERRQ(ierr);
381*c4762a1bSJed Brown       if (key == 1) {
382*c4762a1bSJed Brown 	bus = (VERTEX_Power)(component);
383*c4762a1bSJed Brown 	xarr[offset] = bus->va*PETSC_PI/180.0;
384*c4762a1bSJed Brown 	xarr[offset+1] = bus->vm;
385*c4762a1bSJed Brown       } else if(key == 2) {
386*c4762a1bSJed Brown 	gen = (GEN)(component);
387*c4762a1bSJed Brown 	if (!gen->status) continue;
388*c4762a1bSJed Brown 	xarr[offset+1] = gen->vs;
389*c4762a1bSJed Brown 	break;
390*c4762a1bSJed Brown       }
391*c4762a1bSJed Brown     }
392*c4762a1bSJed Brown   }
393*c4762a1bSJed Brown   ierr = VecRestoreArray(localX,&xarr);CHKERRQ(ierr);
394*c4762a1bSJed Brown   PetscFunctionReturn(0);
395*c4762a1bSJed Brown }
396*c4762a1bSJed Brown 
397*c4762a1bSJed Brown PetscErrorCode SetInitialValues(DM networkdm, Vec X,void* appctx)
398*c4762a1bSJed Brown {
399*c4762a1bSJed Brown   PetscErrorCode ierr;
400*c4762a1bSJed Brown   PetscInt       nv,ne;
401*c4762a1bSJed Brown   const PetscInt *vtx,*edges;
402*c4762a1bSJed Brown   Vec            localX;
403*c4762a1bSJed Brown 
404*c4762a1bSJed Brown   PetscFunctionBegin;
405*c4762a1bSJed Brown   ierr = DMGetLocalVector(networkdm,&localX);CHKERRQ(ierr);
406*c4762a1bSJed Brown 
407*c4762a1bSJed Brown   ierr = VecSet(X,0.0);CHKERRQ(ierr);
408*c4762a1bSJed Brown   ierr = DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
409*c4762a1bSJed Brown   ierr = DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX);CHKERRQ(ierr);
410*c4762a1bSJed Brown 
411*c4762a1bSJed Brown   /* Set initial guess for first subnetwork */
412*c4762a1bSJed Brown   ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
413*c4762a1bSJed Brown   ierr = SetInitialValues_Subnet(networkdm,localX,nv,ne,vtx,edges,appctx);CHKERRQ(ierr);
414*c4762a1bSJed Brown 
415*c4762a1bSJed Brown   /* Set initial guess for second subnetwork */
416*c4762a1bSJed Brown   ierr = DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
417*c4762a1bSJed Brown   ierr = SetInitialValues_Subnet(networkdm,localX,nv,ne,vtx,edges,appctx);CHKERRQ(ierr);
418*c4762a1bSJed Brown 
419*c4762a1bSJed Brown   ierr = DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr);
420*c4762a1bSJed Brown   ierr = DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X);CHKERRQ(ierr);
421*c4762a1bSJed Brown   ierr = DMRestoreLocalVector(networkdm,&localX);CHKERRQ(ierr);
422*c4762a1bSJed Brown   PetscFunctionReturn(0);
423*c4762a1bSJed Brown }
424*c4762a1bSJed Brown 
425*c4762a1bSJed Brown int main(int argc,char ** argv)
426*c4762a1bSJed Brown {
427*c4762a1bSJed Brown   PetscErrorCode   ierr;
428*c4762a1bSJed Brown   char             pfdata_file[PETSC_MAX_PATH_LEN]="case9.m";
429*c4762a1bSJed Brown   PFDATA           *pfdata1,*pfdata2;
430*c4762a1bSJed Brown   PetscInt         numEdges1=0,numVertices1=0,numEdges2=0,numVertices2=0;
431*c4762a1bSJed Brown   PetscInt         *edgelist1 = NULL,*edgelist2 = NULL;
432*c4762a1bSJed Brown   DM               networkdm;
433*c4762a1bSJed Brown   PetscInt         componentkey[4];
434*c4762a1bSJed Brown   UserCtx_Power    User;
435*c4762a1bSJed Brown   PetscLogStage    stage1,stage2;
436*c4762a1bSJed Brown   PetscMPIInt      rank;
437*c4762a1bSJed Brown   PetscInt         nsubnet = 2;
438*c4762a1bSJed Brown   PetscInt         numVertices[2],numEdges[2];
439*c4762a1bSJed Brown   PetscInt         *edgelist[2];
440*c4762a1bSJed Brown   PetscInt         nv,ne;
441*c4762a1bSJed Brown   const PetscInt   *vtx;
442*c4762a1bSJed Brown   const PetscInt   *edges;
443*c4762a1bSJed Brown   PetscInt         i,j;
444*c4762a1bSJed Brown   PetscInt         genj,loadj;
445*c4762a1bSJed Brown   Vec              X,F;
446*c4762a1bSJed Brown   Mat              J;
447*c4762a1bSJed Brown   SNES             snes;
448*c4762a1bSJed Brown 
449*c4762a1bSJed Brown 
450*c4762a1bSJed Brown   ierr = PetscInitialize(&argc,&argv,"poweroptions",help);if (ierr) return ierr;
451*c4762a1bSJed Brown   ierr = MPI_Comm_rank(PETSC_COMM_WORLD,&rank);CHKERRQ(ierr);
452*c4762a1bSJed Brown   {
453*c4762a1bSJed Brown     /* introduce the const crank so the clang static analyzer realizes that if it enters any of the if (crank) then it must have entered the first */
454*c4762a1bSJed Brown     /* this is an experiment to see how the analyzer reacts */
455*c4762a1bSJed Brown     const PetscMPIInt crank = rank;
456*c4762a1bSJed Brown 
457*c4762a1bSJed Brown     /* Create an empty network object */
458*c4762a1bSJed Brown     ierr = DMNetworkCreate(PETSC_COMM_WORLD,&networkdm);CHKERRQ(ierr);
459*c4762a1bSJed Brown 
460*c4762a1bSJed Brown     /* Register the components in the network */
461*c4762a1bSJed Brown     ierr = DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&componentkey[0]);CHKERRQ(ierr);
462*c4762a1bSJed Brown     ierr = DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&componentkey[1]);CHKERRQ(ierr);
463*c4762a1bSJed Brown     ierr = DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&componentkey[2]);CHKERRQ(ierr);
464*c4762a1bSJed Brown     ierr = DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&componentkey[3]);CHKERRQ(ierr);
465*c4762a1bSJed Brown 
466*c4762a1bSJed Brown     ierr = PetscLogStageRegister("Read Data",&stage1);CHKERRQ(ierr);
467*c4762a1bSJed Brown     PetscLogStagePush(stage1);
468*c4762a1bSJed Brown     /* READ THE DATA */
469*c4762a1bSJed Brown     if (!crank) {
470*c4762a1bSJed Brown       /* Only rank 0 reads the data */
471*c4762a1bSJed Brown       ierr = PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,PETSC_MAX_PATH_LEN-1,NULL);CHKERRQ(ierr);
472*c4762a1bSJed Brown       /* HERE WE CREATE COPIES OF THE SAME NETWORK THAT WILL BE TREATED AS SUBNETWORKS */
473*c4762a1bSJed Brown 
474*c4762a1bSJed Brown       /*    READ DATA FOR THE FIRST SUBNETWORK */
475*c4762a1bSJed Brown       ierr = PetscNew(&pfdata1);CHKERRQ(ierr);
476*c4762a1bSJed Brown       ierr = PFReadMatPowerData(pfdata1,pfdata_file);CHKERRQ(ierr);
477*c4762a1bSJed Brown       User.Sbase = pfdata1->sbase;
478*c4762a1bSJed Brown 
479*c4762a1bSJed Brown       numEdges1 = pfdata1->nbranch;
480*c4762a1bSJed Brown       numVertices1 = pfdata1->nbus;
481*c4762a1bSJed Brown 
482*c4762a1bSJed Brown       ierr = PetscMalloc1(2*numEdges1,&edgelist1);CHKERRQ(ierr);
483*c4762a1bSJed Brown       ierr = GetListofEdges_Power(pfdata1,edgelist1);CHKERRQ(ierr);
484*c4762a1bSJed Brown 
485*c4762a1bSJed Brown       /*    READ DATA FOR THE SECOND SUBNETWORK */
486*c4762a1bSJed Brown       ierr = PetscNew(&pfdata2);CHKERRQ(ierr);
487*c4762a1bSJed Brown       ierr = PFReadMatPowerData(pfdata2,pfdata_file);CHKERRQ(ierr);
488*c4762a1bSJed Brown       User.Sbase = pfdata2->sbase;
489*c4762a1bSJed Brown 
490*c4762a1bSJed Brown       numEdges2 = pfdata2->nbranch;
491*c4762a1bSJed Brown       numVertices2 = pfdata2->nbus;
492*c4762a1bSJed Brown 
493*c4762a1bSJed Brown       ierr = PetscMalloc1(2*numEdges2,&edgelist2);CHKERRQ(ierr);
494*c4762a1bSJed Brown       ierr = GetListofEdges_Power(pfdata2,edgelist2);CHKERRQ(ierr);
495*c4762a1bSJed Brown     }
496*c4762a1bSJed Brown 
497*c4762a1bSJed Brown     PetscLogStagePop();
498*c4762a1bSJed Brown     ierr = MPI_Barrier(PETSC_COMM_WORLD);CHKERRQ(ierr);
499*c4762a1bSJed Brown     ierr = PetscLogStageRegister("Create network",&stage2);CHKERRQ(ierr);
500*c4762a1bSJed Brown     PetscLogStagePush(stage2);
501*c4762a1bSJed Brown 
502*c4762a1bSJed Brown     /* Set number of nodes/edges */
503*c4762a1bSJed Brown     numVertices[0] = numVertices1; numVertices[1] = numVertices2;
504*c4762a1bSJed Brown     numEdges[0] = numEdges1; numEdges[1] = numEdges2;
505*c4762a1bSJed Brown     ierr = DMNetworkSetSizes(networkdm,nsubnet,numVertices,numEdges,0,NULL);CHKERRQ(ierr);
506*c4762a1bSJed Brown 
507*c4762a1bSJed Brown     edgelist[0] = edgelist1; edgelist[1] = edgelist2;
508*c4762a1bSJed Brown 
509*c4762a1bSJed Brown     /* Add edge connectivity */
510*c4762a1bSJed Brown     ierr = DMNetworkSetEdgeList(networkdm,edgelist,NULL);CHKERRQ(ierr);
511*c4762a1bSJed Brown 
512*c4762a1bSJed Brown     /* Set up the network layout */
513*c4762a1bSJed Brown     ierr = DMNetworkLayoutSetUp(networkdm);CHKERRQ(ierr);
514*c4762a1bSJed Brown 
515*c4762a1bSJed Brown     /* Add network components only process 0 has any data to add*/
516*c4762a1bSJed Brown     if (!crank) {
517*c4762a1bSJed Brown       genj=0; loadj=0;
518*c4762a1bSJed Brown 
519*c4762a1bSJed Brown       /* ADD VARIABLES AND COMPONENTS FOR THE FIRST SUBNETWORK */
520*c4762a1bSJed Brown       ierr = DMNetworkGetSubnetworkInfo(networkdm,0,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
521*c4762a1bSJed Brown 
522*c4762a1bSJed Brown       for (i = 0; i < ne; i++) {
523*c4762a1bSJed Brown         ierr = DMNetworkAddComponent(networkdm,edges[i],componentkey[0],&pfdata1->branch[i]);CHKERRQ(ierr);
524*c4762a1bSJed Brown       }
525*c4762a1bSJed Brown 
526*c4762a1bSJed Brown       for (i = 0; i < nv; i++) {
527*c4762a1bSJed Brown         ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[1],&pfdata1->bus[i]);CHKERRQ(ierr);
528*c4762a1bSJed Brown         if (pfdata1->bus[i].ngen) {
529*c4762a1bSJed Brown           for (j = 0; j < pfdata1->bus[i].ngen; j++) {
530*c4762a1bSJed Brown             ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[2],&pfdata1->gen[genj++]);CHKERRQ(ierr);
531*c4762a1bSJed Brown           }
532*c4762a1bSJed Brown         }
533*c4762a1bSJed Brown         if (pfdata1->bus[i].nload) {
534*c4762a1bSJed Brown           for (j=0; j < pfdata1->bus[i].nload; j++) {
535*c4762a1bSJed Brown             ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[3],&pfdata1->load[loadj++]);CHKERRQ(ierr);
536*c4762a1bSJed Brown           }
537*c4762a1bSJed Brown         }
538*c4762a1bSJed Brown         /* Add number of variables */
539*c4762a1bSJed Brown         ierr = DMNetworkAddNumVariables(networkdm,vtx[i],2);CHKERRQ(ierr);
540*c4762a1bSJed Brown       }
541*c4762a1bSJed Brown 
542*c4762a1bSJed Brown       genj=0; loadj=0;
543*c4762a1bSJed Brown 
544*c4762a1bSJed Brown       /* ADD VARIABLES AND COMPONENTS FOR THE SECOND SUBNETWORK */
545*c4762a1bSJed Brown       ierr = DMNetworkGetSubnetworkInfo(networkdm,1,&nv,&ne,&vtx,&edges);CHKERRQ(ierr);
546*c4762a1bSJed Brown 
547*c4762a1bSJed Brown       for (i = 0; i < ne; i++) {
548*c4762a1bSJed Brown         ierr = DMNetworkAddComponent(networkdm,edges[i],componentkey[0],&pfdata2->branch[i]);CHKERRQ(ierr);
549*c4762a1bSJed Brown       }
550*c4762a1bSJed Brown 
551*c4762a1bSJed Brown       for (i = 0; i < nv; i++) {
552*c4762a1bSJed Brown         ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[1],&pfdata2->bus[i]);CHKERRQ(ierr);
553*c4762a1bSJed Brown         if (pfdata2->bus[i].ngen) {
554*c4762a1bSJed Brown           for (j = 0; j < pfdata2->bus[i].ngen; j++) {
555*c4762a1bSJed Brown             ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[2],&pfdata2->gen[genj++]);CHKERRQ(ierr);
556*c4762a1bSJed Brown           }
557*c4762a1bSJed Brown         }
558*c4762a1bSJed Brown         if (pfdata2->bus[i].nload) {
559*c4762a1bSJed Brown           for (j=0; j < pfdata2->bus[i].nload; j++) {
560*c4762a1bSJed Brown             ierr = DMNetworkAddComponent(networkdm,vtx[i],componentkey[3],&pfdata2->load[loadj++]);CHKERRQ(ierr);
561*c4762a1bSJed Brown           }
562*c4762a1bSJed Brown         }
563*c4762a1bSJed Brown         /* Add number of variables */
564*c4762a1bSJed Brown         ierr = DMNetworkAddNumVariables(networkdm,vtx[i],2);CHKERRQ(ierr);
565*c4762a1bSJed Brown       }
566*c4762a1bSJed Brown     }
567*c4762a1bSJed Brown 
568*c4762a1bSJed Brown     /* Set up DM for use */
569*c4762a1bSJed Brown     ierr = DMSetUp(networkdm);CHKERRQ(ierr);
570*c4762a1bSJed Brown 
571*c4762a1bSJed Brown     if (!crank) {
572*c4762a1bSJed Brown       ierr = PetscFree(edgelist1);CHKERRQ(ierr);
573*c4762a1bSJed Brown       ierr = PetscFree(edgelist2);CHKERRQ(ierr);
574*c4762a1bSJed Brown     }
575*c4762a1bSJed Brown 
576*c4762a1bSJed Brown     if (!crank) {
577*c4762a1bSJed Brown       ierr = PetscFree(pfdata1->bus);CHKERRQ(ierr);
578*c4762a1bSJed Brown       ierr = PetscFree(pfdata1->gen);CHKERRQ(ierr);
579*c4762a1bSJed Brown       ierr = PetscFree(pfdata1->branch);CHKERRQ(ierr);
580*c4762a1bSJed Brown       ierr = PetscFree(pfdata1->load);CHKERRQ(ierr);
581*c4762a1bSJed Brown       ierr = PetscFree(pfdata1);CHKERRQ(ierr);
582*c4762a1bSJed Brown 
583*c4762a1bSJed Brown       ierr = PetscFree(pfdata2->bus);CHKERRQ(ierr);
584*c4762a1bSJed Brown       ierr = PetscFree(pfdata2->gen);CHKERRQ(ierr);
585*c4762a1bSJed Brown       ierr = PetscFree(pfdata2->branch);CHKERRQ(ierr);
586*c4762a1bSJed Brown       ierr = PetscFree(pfdata2->load);CHKERRQ(ierr);
587*c4762a1bSJed Brown       ierr = PetscFree(pfdata2);CHKERRQ(ierr);
588*c4762a1bSJed Brown     }
589*c4762a1bSJed Brown 
590*c4762a1bSJed Brown     /* Distribute networkdm to multiple processes */
591*c4762a1bSJed Brown     ierr = DMNetworkDistribute(&networkdm,0);CHKERRQ(ierr);
592*c4762a1bSJed Brown 
593*c4762a1bSJed Brown     PetscLogStagePop();
594*c4762a1bSJed Brown 
595*c4762a1bSJed Brown     /* Broadcast Sbase to all processors */
596*c4762a1bSJed Brown     ierr = MPI_Bcast(&User.Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD);CHKERRQ(ierr);
597*c4762a1bSJed Brown 
598*c4762a1bSJed Brown     ierr = DMCreateGlobalVector(networkdm,&X);CHKERRQ(ierr);
599*c4762a1bSJed Brown     ierr = VecDuplicate(X,&F);CHKERRQ(ierr);
600*c4762a1bSJed Brown 
601*c4762a1bSJed Brown     ierr = DMCreateMatrix(networkdm,&J);CHKERRQ(ierr);
602*c4762a1bSJed Brown     ierr = MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE);CHKERRQ(ierr);
603*c4762a1bSJed Brown 
604*c4762a1bSJed Brown     ierr = SetInitialValues(networkdm,X,&User);CHKERRQ(ierr);
605*c4762a1bSJed Brown 
606*c4762a1bSJed Brown     /* HOOK UP SOLVER */
607*c4762a1bSJed Brown     ierr = SNESCreate(PETSC_COMM_WORLD,&snes);CHKERRQ(ierr);
608*c4762a1bSJed Brown     ierr = SNESSetDM(snes,networkdm);CHKERRQ(ierr);
609*c4762a1bSJed Brown     ierr = SNESSetFunction(snes,F,FormFunction,&User);CHKERRQ(ierr);
610*c4762a1bSJed Brown     ierr = SNESSetJacobian(snes,J,J,FormJacobian,&User);CHKERRQ(ierr);
611*c4762a1bSJed Brown     ierr = SNESSetFromOptions(snes);CHKERRQ(ierr);
612*c4762a1bSJed Brown 
613*c4762a1bSJed Brown     ierr = SNESSolve(snes,NULL,X);CHKERRQ(ierr);
614*c4762a1bSJed Brown     /* ierr = VecView(X,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); */
615*c4762a1bSJed Brown 
616*c4762a1bSJed Brown     ierr = VecDestroy(&X);CHKERRQ(ierr);
617*c4762a1bSJed Brown     ierr = VecDestroy(&F);CHKERRQ(ierr);
618*c4762a1bSJed Brown     ierr = MatDestroy(&J);CHKERRQ(ierr);
619*c4762a1bSJed Brown 
620*c4762a1bSJed Brown     ierr = SNESDestroy(&snes);CHKERRQ(ierr);
621*c4762a1bSJed Brown     ierr = DMDestroy(&networkdm);CHKERRQ(ierr);
622*c4762a1bSJed Brown   }
623*c4762a1bSJed Brown   ierr = PetscFinalize();
624*c4762a1bSJed Brown   return ierr;
625*c4762a1bSJed Brown }
626*c4762a1bSJed Brown 
627*c4762a1bSJed Brown /*TEST
628*c4762a1bSJed Brown 
629*c4762a1bSJed Brown    build:
630*c4762a1bSJed Brown      depends: PFReadData.c pffunctions.c
631*c4762a1bSJed Brown      requires: !complex double define(PETSC_HAVE_ATTRIBUTEALIGNED)
632*c4762a1bSJed Brown 
633*c4762a1bSJed Brown 
634*c4762a1bSJed Brown    test:
635*c4762a1bSJed Brown      args: -snes_rtol 1.e-3
636*c4762a1bSJed Brown      localrunfiles: poweroptions case9.m
637*c4762a1bSJed Brown      output_file: output/power2_1.out
638*c4762a1bSJed Brown      requires: double !complex
639*c4762a1bSJed Brown 
640*c4762a1bSJed Brown    test:
641*c4762a1bSJed Brown      suffix: 2
642*c4762a1bSJed Brown      args: -snes_rtol 1.e-3 -petscpartitioner_type simple
643*c4762a1bSJed Brown      nsize: 4
644*c4762a1bSJed Brown      localrunfiles: poweroptions case9.m
645*c4762a1bSJed Brown      output_file: output/power2_1.out
646*c4762a1bSJed Brown      requires: double !complex
647*c4762a1bSJed Brown 
648*c4762a1bSJed Brown TEST*/
649