xref: /petsc/src/snes/tutorials/network/power/power2.c (revision 327415f76d85372a4417cf1aaa14db707d4d6c04)
1c4762a1bSJed Brown static char help[] = "This example demonstrates the use of DMNetwork interface with subnetworks for solving a nonlinear electric power grid problem.\n\
2c4762a1bSJed Brown                       The available solver options are in the poweroptions file and the data files are in the datafiles directory.\n\
3c4762a1bSJed Brown                       The data file format used is from the MatPower package (http://www.pserc.cornell.edu//matpower/).\n\
4c4762a1bSJed Brown                       This example shows the use of subnetwork feature in DMNetwork. It creates duplicates of the same network which are treated as subnetworks.\n\
5c4762a1bSJed Brown                       Run this program: mpiexec -n <n> ./pf2\n\
6c4762a1bSJed Brown                       mpiexec -n <n> ./pf2 \n";
7c4762a1bSJed Brown 
8c4762a1bSJed Brown #include "power.h"
9c4762a1bSJed Brown #include <petscdmnetwork.h>
10c4762a1bSJed Brown 
11c4762a1bSJed Brown PetscErrorCode FormFunction_Subnet(DM networkdm,Vec localX, Vec localF,PetscInt nv,PetscInt ne,const PetscInt* vtx,const PetscInt* edges,void* appctx)
12c4762a1bSJed Brown {
13c4762a1bSJed Brown   UserCtx_Power     *User = (UserCtx_Power*)appctx;
14c4762a1bSJed Brown   PetscInt          e,v,vfrom,vto;
15c4762a1bSJed Brown   const PetscScalar *xarr;
16c4762a1bSJed Brown   PetscScalar       *farr;
17c4762a1bSJed Brown   PetscInt          offsetfrom,offsetto,offset;
18c4762a1bSJed Brown 
19c4762a1bSJed Brown   PetscFunctionBegin;
209566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localX,&xarr));
219566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localF,&farr));
22c4762a1bSJed Brown 
23c4762a1bSJed Brown   for (v=0; v<nv; v++) {
24c4762a1bSJed Brown     PetscInt      i,j,key;
25c4762a1bSJed Brown     PetscScalar   Vm;
26c4762a1bSJed Brown     PetscScalar   Sbase = User->Sbase;
27c4762a1bSJed Brown     VERTEX_Power  bus = NULL;
28c4762a1bSJed Brown     GEN           gen;
29c4762a1bSJed Brown     LOAD          load;
30c4762a1bSJed Brown     PetscBool     ghostvtex;
31c4762a1bSJed Brown     PetscInt      numComps;
32c4762a1bSJed Brown     void*         component;
33c4762a1bSJed Brown 
349566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex));
359566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm,vtx[v],&numComps));
369566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset));
37c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
389566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL));
39c4762a1bSJed Brown       if (key == 1) {
40c4762a1bSJed Brown         PetscInt       nconnedges;
41c4762a1bSJed Brown         const PetscInt *connedges;
42c4762a1bSJed Brown 
43c4762a1bSJed Brown         bus = (VERTEX_Power)(component);
44c4762a1bSJed Brown         /* Handle reference bus constrained dofs */
45c4762a1bSJed Brown         if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
46c4762a1bSJed Brown           farr[offset] = xarr[offset] - bus->va*PETSC_PI/180.0;
47c4762a1bSJed Brown           farr[offset+1] = xarr[offset+1] - bus->vm;
48c4762a1bSJed Brown           break;
49c4762a1bSJed Brown         }
50c4762a1bSJed Brown 
51c4762a1bSJed Brown         if (!ghostvtex) {
52c4762a1bSJed Brown           Vm = xarr[offset+1];
53c4762a1bSJed Brown 
54c4762a1bSJed Brown           /* Shunt injections */
55c4762a1bSJed Brown           farr[offset] += Vm*Vm*bus->gl/Sbase;
56c4762a1bSJed Brown           if (bus->ide != PV_BUS) farr[offset+1] += -Vm*Vm*bus->bl/Sbase;
57c4762a1bSJed Brown         }
58c4762a1bSJed Brown 
599566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges));
60c4762a1bSJed Brown         for (i=0; i < nconnedges; i++) {
61c4762a1bSJed Brown           EDGE_Power     branch;
62c4762a1bSJed Brown           PetscInt       keye;
63c4762a1bSJed Brown           PetscScalar    Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt;
64c4762a1bSJed Brown           const PetscInt *cone;
65c4762a1bSJed Brown           PetscScalar    Vmf,Vmt,thetaf,thetat,thetaft,thetatf;
66c4762a1bSJed Brown 
67c4762a1bSJed Brown           e = connedges[i];
689566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetComponent(networkdm,e,0,&keye,(void**)&branch,NULL));
69c4762a1bSJed Brown           if (!branch->status) continue;
70c4762a1bSJed Brown           Gff = branch->yff[0];
71c4762a1bSJed Brown           Bff = branch->yff[1];
72c4762a1bSJed Brown           Gft = branch->yft[0];
73c4762a1bSJed Brown           Bft = branch->yft[1];
74c4762a1bSJed Brown           Gtf = branch->ytf[0];
75c4762a1bSJed Brown           Btf = branch->ytf[1];
76c4762a1bSJed Brown           Gtt = branch->ytt[0];
77c4762a1bSJed Brown           Btt = branch->ytt[1];
78c4762a1bSJed Brown 
799566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetConnectedVertices(networkdm,e,&cone));
80c4762a1bSJed Brown           vfrom = cone[0];
81c4762a1bSJed Brown           vto   = cone[1];
82c4762a1bSJed Brown 
839566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom));
849566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto));
85c4762a1bSJed Brown 
86c4762a1bSJed Brown           thetaf = xarr[offsetfrom];
87c4762a1bSJed Brown           Vmf     = xarr[offsetfrom+1];
88c4762a1bSJed Brown           thetat  = xarr[offsetto];
89c4762a1bSJed Brown           Vmt     = xarr[offsetto+1];
90c4762a1bSJed Brown           thetaft = thetaf - thetat;
91c4762a1bSJed Brown           thetatf = thetat - thetaf;
92c4762a1bSJed Brown 
93c4762a1bSJed Brown           if (vfrom == vtx[v]) {
94c4762a1bSJed Brown             farr[offsetfrom]   += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft));
95c4762a1bSJed Brown             farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
96c4762a1bSJed Brown           } else {
97c4762a1bSJed Brown             farr[offsetto]   += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf));
98c4762a1bSJed Brown             farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
99c4762a1bSJed Brown           }
100c4762a1bSJed Brown         }
101c4762a1bSJed Brown       } else if (key == 2) {
102c4762a1bSJed Brown         if (!ghostvtex) {
103c4762a1bSJed Brown           gen = (GEN)(component);
104c4762a1bSJed Brown           if (!gen->status) continue;
105c4762a1bSJed Brown           farr[offset] += -gen->pg/Sbase;
106c4762a1bSJed Brown           farr[offset+1] += -gen->qg/Sbase;
107c4762a1bSJed Brown         }
108c4762a1bSJed Brown       } else if (key == 3) {
109c4762a1bSJed Brown         if (!ghostvtex) {
110c4762a1bSJed Brown           load = (LOAD)(component);
111c4762a1bSJed Brown           farr[offset] += load->pl/Sbase;
112c4762a1bSJed Brown           farr[offset+1] += load->ql/Sbase;
113c4762a1bSJed Brown         }
114c4762a1bSJed Brown       }
115c4762a1bSJed Brown     }
116c4762a1bSJed Brown     if (bus && bus->ide == PV_BUS) {
117c4762a1bSJed Brown       farr[offset+1] = xarr[offset+1] - bus->vm;
118c4762a1bSJed Brown     }
119c4762a1bSJed Brown   }
1209566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localX,&xarr));
1219566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localF,&farr));
122c4762a1bSJed Brown   PetscFunctionReturn(0);
123c4762a1bSJed Brown }
124c4762a1bSJed Brown 
125c4762a1bSJed Brown PetscErrorCode FormFunction(SNES snes,Vec X, Vec F,void *appctx)
126c4762a1bSJed Brown {
127c4762a1bSJed Brown   DM             networkdm;
128c4762a1bSJed Brown   Vec            localX,localF;
129c4762a1bSJed Brown   PetscInt       nv,ne;
130c4762a1bSJed Brown   const PetscInt *vtx,*edges;
131c4762a1bSJed Brown 
132c4762a1bSJed Brown   PetscFunctionBegin;
1339566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&networkdm));
1349566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm,&localX));
1359566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm,&localF));
1369566063dSJacob Faibussowitsch   PetscCall(VecSet(F,0.0));
137c4762a1bSJed Brown 
1389566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX));
1399566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX));
140c4762a1bSJed Brown 
1419566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm,F,INSERT_VALUES,localF));
1429566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm,F,INSERT_VALUES,localF));
143c4762a1bSJed Brown 
144c4762a1bSJed Brown   /* Form Function for first subnetwork */
1459566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges));
1469566063dSJacob Faibussowitsch   PetscCall(FormFunction_Subnet(networkdm,localX,localF,nv,ne,vtx,edges,appctx));
147c4762a1bSJed Brown 
148c4762a1bSJed Brown   /* Form Function for second subnetwork */
1499566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges));
1509566063dSJacob Faibussowitsch   PetscCall(FormFunction_Subnet(networkdm,localX,localF,nv,ne,vtx,edges,appctx));
151c4762a1bSJed Brown 
1529566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm,&localX));
153c4762a1bSJed Brown 
1549566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(networkdm,localF,ADD_VALUES,F));
1559566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(networkdm,localF,ADD_VALUES,F));
1569566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm,&localF));
157c4762a1bSJed Brown   PetscFunctionReturn(0);
158c4762a1bSJed Brown }
159c4762a1bSJed Brown 
160c4762a1bSJed Brown PetscErrorCode FormJacobian_Subnet(DM networkdm,Vec localX, Mat J, Mat Jpre, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
161c4762a1bSJed Brown {
162c4762a1bSJed Brown   UserCtx_Power     *User=(UserCtx_Power*)appctx;
163c4762a1bSJed Brown   PetscInt          e,v,vfrom,vto;
164c4762a1bSJed Brown   const PetscScalar *xarr;
165c4762a1bSJed Brown   PetscInt          offsetfrom,offsetto,goffsetfrom,goffsetto;
166c4762a1bSJed Brown   PetscInt          row[2],col[8];
167c4762a1bSJed Brown   PetscScalar       values[8];
168c4762a1bSJed Brown 
169c4762a1bSJed Brown   PetscFunctionBegin;
1709566063dSJacob Faibussowitsch   PetscCall(VecGetArrayRead(localX,&xarr));
171c4762a1bSJed Brown 
172c4762a1bSJed Brown   for (v=0; v<nv; v++) {
173c4762a1bSJed Brown     PetscInt    i,j,key;
174c4762a1bSJed Brown     PetscInt    offset,goffset;
175c4762a1bSJed Brown     PetscScalar Vm;
176c4762a1bSJed Brown     PetscScalar Sbase=User->Sbase;
177c4762a1bSJed Brown     VERTEX_Power bus;
178c4762a1bSJed Brown     PetscBool   ghostvtex;
179c4762a1bSJed Brown     PetscInt    numComps;
180c4762a1bSJed Brown     void*       component;
181c4762a1bSJed Brown 
1829566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm,vtx[v],&ghostvtex));
1839566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm,vtx[v],&numComps));
184c4762a1bSJed Brown     for (j = 0; j < numComps; j++) {
1859566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetLocalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&offset));
1869566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetGlobalVecOffset(networkdm,vtx[v],ALL_COMPONENTS,&goffset));
1879566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(networkdm,vtx[v],j,&key,&component,NULL));
188c4762a1bSJed Brown       if (key == 1) {
189c4762a1bSJed Brown         PetscInt       nconnedges;
190c4762a1bSJed Brown         const PetscInt *connedges;
191c4762a1bSJed Brown 
192c4762a1bSJed Brown         bus = (VERTEX_Power)(component);
193c4762a1bSJed Brown         if (!ghostvtex) {
194c4762a1bSJed Brown           /* Handle reference bus constrained dofs */
195c4762a1bSJed Brown           if (bus->ide == REF_BUS || bus->ide == ISOLATED_BUS) {
196c4762a1bSJed Brown             row[0] = goffset; row[1] = goffset+1;
197c4762a1bSJed Brown             col[0] = goffset; col[1] = goffset+1; col[2] = goffset; col[3] = goffset+1;
198c4762a1bSJed Brown             values[0] = 1.0; values[1] = 0.0; values[2] = 0.0; values[3] = 1.0;
1999566063dSJacob Faibussowitsch             PetscCall(MatSetValues(J,2,row,2,col,values,ADD_VALUES));
200c4762a1bSJed Brown             break;
201c4762a1bSJed Brown           }
202c4762a1bSJed Brown 
203c4762a1bSJed Brown           Vm = xarr[offset+1];
204c4762a1bSJed Brown 
205c4762a1bSJed Brown           /* Shunt injections */
206c4762a1bSJed Brown           row[0] = goffset; row[1] = goffset+1;
207c4762a1bSJed Brown           col[0] = goffset; col[1] = goffset+1;
208c4762a1bSJed Brown           values[0] = values[1] = values[2] = values[3] = 0.0;
209c4762a1bSJed Brown           if (bus->ide != PV_BUS) {
210c4762a1bSJed Brown             values[1] = 2.0*Vm*bus->gl/Sbase;
211c4762a1bSJed Brown             values[3] = -2.0*Vm*bus->bl/Sbase;
212c4762a1bSJed Brown           }
2139566063dSJacob Faibussowitsch           PetscCall(MatSetValues(J,2,row,2,col,values,ADD_VALUES));
214c4762a1bSJed Brown         }
215c4762a1bSJed Brown 
2169566063dSJacob Faibussowitsch         PetscCall(DMNetworkGetSupportingEdges(networkdm,vtx[v],&nconnedges,&connedges));
217c4762a1bSJed Brown         for (i=0; i < nconnedges; i++) {
218c4762a1bSJed Brown           EDGE_Power       branch;
219c4762a1bSJed Brown           VERTEX_Power     busf,bust;
220c4762a1bSJed Brown           PetscInt       keyf,keyt;
221c4762a1bSJed Brown           PetscScalar    Gff,Bff,Gft,Bft,Gtf,Btf,Gtt,Btt;
222c4762a1bSJed Brown           const PetscInt *cone;
223c4762a1bSJed Brown           PetscScalar    Vmf,Vmt,thetaf,thetat,thetaft,thetatf;
224c4762a1bSJed Brown 
225c4762a1bSJed Brown           e = connedges[i];
2269566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetComponent(networkdm,e,0,&key,(void**)&branch,NULL));
227c4762a1bSJed Brown           if (!branch->status) continue;
228c4762a1bSJed Brown 
229c4762a1bSJed Brown           Gff = branch->yff[0];
230c4762a1bSJed Brown           Bff = branch->yff[1];
231c4762a1bSJed Brown           Gft = branch->yft[0];
232c4762a1bSJed Brown           Bft = branch->yft[1];
233c4762a1bSJed Brown           Gtf = branch->ytf[0];
234c4762a1bSJed Brown           Btf = branch->ytf[1];
235c4762a1bSJed Brown           Gtt = branch->ytt[0];
236c4762a1bSJed Brown           Btt = branch->ytt[1];
237c4762a1bSJed Brown 
2389566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetConnectedVertices(networkdm,e,&cone));
239c4762a1bSJed Brown           vfrom = cone[0];
240c4762a1bSJed Brown           vto   = cone[1];
241c4762a1bSJed Brown 
2429566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&offsetfrom));
2439566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetLocalVecOffset(networkdm,vto,ALL_COMPONENTS,&offsetto));
2449566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetGlobalVecOffset(networkdm,vfrom,ALL_COMPONENTS,&goffsetfrom));
2459566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetGlobalVecOffset(networkdm,vto,ALL_COMPONENTS,&goffsetto));
246c4762a1bSJed Brown 
247c4762a1bSJed Brown           if (goffsetto < 0) goffsetto = -goffsetto - 1;
248c4762a1bSJed Brown 
249c4762a1bSJed Brown           thetaf = xarr[offsetfrom];
250c4762a1bSJed Brown           Vmf     = xarr[offsetfrom+1];
251c4762a1bSJed Brown           thetat = xarr[offsetto];
252c4762a1bSJed Brown           Vmt     = xarr[offsetto+1];
253c4762a1bSJed Brown           thetaft = thetaf - thetat;
254c4762a1bSJed Brown           thetatf = thetat - thetaf;
255c4762a1bSJed Brown 
2569566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetComponent(networkdm,vfrom,0,&keyf,(void**)&busf,NULL));
2579566063dSJacob Faibussowitsch           PetscCall(DMNetworkGetComponent(networkdm,vto,0,&keyt,(void**)&bust,NULL));
258c4762a1bSJed Brown 
259c4762a1bSJed Brown           if (vfrom == vtx[v]) {
260c4762a1bSJed Brown             if (busf->ide != REF_BUS) {
261c4762a1bSJed Brown               /*    farr[offsetfrom]   += Gff*Vmf*Vmf + Vmf*Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft));  */
262c4762a1bSJed Brown               row[0]  = goffsetfrom;
263c4762a1bSJed Brown               col[0]  = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1;
264c4762a1bSJed Brown               values[0] =  Vmf*Vmt*(Gft*-PetscSinScalar(thetaft) + Bft*PetscCosScalar(thetaft)); /* df_dthetaf */
265c4762a1bSJed Brown               values[1] =  2.0*Gff*Vmf + Vmt*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmf */
266c4762a1bSJed Brown               values[2] =  Vmf*Vmt*(Gft*PetscSinScalar(thetaft) + Bft*-PetscCosScalar(thetaft)); /* df_dthetat */
267c4762a1bSJed Brown               values[3] =  Vmf*(Gft*PetscCosScalar(thetaft) + Bft*PetscSinScalar(thetaft)); /* df_dVmt */
268c4762a1bSJed Brown 
2699566063dSJacob Faibussowitsch               PetscCall(MatSetValues(J,1,row,4,col,values,ADD_VALUES));
270c4762a1bSJed Brown             }
271c4762a1bSJed Brown             if (busf->ide != PV_BUS && busf->ide != REF_BUS) {
272c4762a1bSJed Brown               row[0] = goffsetfrom+1;
273c4762a1bSJed Brown               col[0]  = goffsetfrom; col[1] = goffsetfrom+1; col[2] = goffsetto; col[3] = goffsetto+1;
274c4762a1bSJed Brown               /*    farr[offsetfrom+1] += -Bff*Vmf*Vmf + Vmf*Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft)); */
275c4762a1bSJed Brown               values[0] =  Vmf*Vmt*(Bft*PetscSinScalar(thetaft) + Gft*PetscCosScalar(thetaft));
276c4762a1bSJed Brown               values[1] =  -2.0*Bff*Vmf + Vmt*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
277c4762a1bSJed Brown               values[2] =  Vmf*Vmt*(-Bft*PetscSinScalar(thetaft) + Gft*-PetscCosScalar(thetaft));
278c4762a1bSJed Brown               values[3] =  Vmf*(-Bft*PetscCosScalar(thetaft) + Gft*PetscSinScalar(thetaft));
279c4762a1bSJed Brown 
2809566063dSJacob Faibussowitsch               PetscCall(MatSetValues(J,1,row,4,col,values,ADD_VALUES));
281c4762a1bSJed Brown             }
282c4762a1bSJed Brown           } else {
283c4762a1bSJed Brown             if (bust->ide != REF_BUS) {
284c4762a1bSJed Brown               row[0] = goffsetto;
285c4762a1bSJed Brown               col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1;
286c4762a1bSJed Brown               /*    farr[offsetto]   += Gtt*Vmt*Vmt + Vmt*Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); */
287c4762a1bSJed Brown               values[0] =  Vmt*Vmf*(Gtf*-PetscSinScalar(thetatf) + Btf*PetscCosScalar(thetaft)); /* df_dthetat */
288c4762a1bSJed Brown               values[1] =  2.0*Gtt*Vmt + Vmf*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmt */
289c4762a1bSJed Brown               values[2] =  Vmt*Vmf*(Gtf*PetscSinScalar(thetatf) + Btf*-PetscCosScalar(thetatf)); /* df_dthetaf */
290c4762a1bSJed Brown               values[3] =  Vmt*(Gtf*PetscCosScalar(thetatf) + Btf*PetscSinScalar(thetatf)); /* df_dVmf */
291c4762a1bSJed Brown 
2929566063dSJacob Faibussowitsch               PetscCall(MatSetValues(J,1,row,4,col,values,ADD_VALUES));
293c4762a1bSJed Brown             }
294c4762a1bSJed Brown             if (bust->ide != PV_BUS && bust->ide != REF_BUS) {
295c4762a1bSJed Brown               row[0] = goffsetto+1;
296c4762a1bSJed Brown               col[0] = goffsetto; col[1] = goffsetto+1; col[2] = goffsetfrom; col[3] = goffsetfrom+1;
297c4762a1bSJed Brown               /*    farr[offsetto+1] += -Btt*Vmt*Vmt + Vmt*Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf)); */
298c4762a1bSJed Brown               values[0] =  Vmt*Vmf*(Btf*PetscSinScalar(thetatf) + Gtf*PetscCosScalar(thetatf));
299c4762a1bSJed Brown               values[1] =  -2.0*Btt*Vmt + Vmf*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
300c4762a1bSJed Brown               values[2] =  Vmt*Vmf*(-Btf*PetscSinScalar(thetatf) + Gtf*-PetscCosScalar(thetatf));
301c4762a1bSJed Brown               values[3] =  Vmt*(-Btf*PetscCosScalar(thetatf) + Gtf*PetscSinScalar(thetatf));
302c4762a1bSJed Brown 
3039566063dSJacob Faibussowitsch               PetscCall(MatSetValues(J,1,row,4,col,values,ADD_VALUES));
304c4762a1bSJed Brown             }
305c4762a1bSJed Brown           }
306c4762a1bSJed Brown         }
307c4762a1bSJed Brown         if (!ghostvtex && bus->ide == PV_BUS) {
308c4762a1bSJed Brown           row[0] = goffset+1; col[0] = goffset+1;
309c4762a1bSJed Brown           values[0]  = 1.0;
3109566063dSJacob Faibussowitsch           PetscCall(MatSetValues(J,1,row,1,col,values,ADD_VALUES));
311c4762a1bSJed Brown         }
312c4762a1bSJed Brown       }
313c4762a1bSJed Brown     }
314c4762a1bSJed Brown   }
3159566063dSJacob Faibussowitsch   PetscCall(VecRestoreArrayRead(localX,&xarr));
316c4762a1bSJed Brown   PetscFunctionReturn(0);
317c4762a1bSJed Brown }
318c4762a1bSJed Brown 
319c4762a1bSJed Brown PetscErrorCode FormJacobian(SNES snes,Vec X, Mat J,Mat Jpre,void *appctx)
320c4762a1bSJed Brown {
321c4762a1bSJed Brown   DM             networkdm;
322c4762a1bSJed Brown   Vec            localX;
323c4762a1bSJed Brown   PetscInt       ne,nv;
324c4762a1bSJed Brown   const PetscInt *vtx,*edges;
325c4762a1bSJed Brown 
326c4762a1bSJed Brown   PetscFunctionBegin;
3279566063dSJacob Faibussowitsch   PetscCall(MatZeroEntries(J));
328c4762a1bSJed Brown 
3299566063dSJacob Faibussowitsch   PetscCall(SNESGetDM(snes,&networkdm));
3309566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm,&localX));
331c4762a1bSJed Brown 
3329566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX));
3339566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX));
334c4762a1bSJed Brown 
335c4762a1bSJed Brown   /* Form Jacobian for first subnetwork */
3369566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges));
3379566063dSJacob Faibussowitsch   PetscCall(FormJacobian_Subnet(networkdm,localX,J,Jpre,nv,ne,vtx,edges,appctx));
338c4762a1bSJed Brown 
339c4762a1bSJed Brown   /* Form Jacobian for second subnetwork */
3409566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges));
3419566063dSJacob Faibussowitsch   PetscCall(FormJacobian_Subnet(networkdm,localX,J,Jpre,nv,ne,vtx,edges,appctx));
342c4762a1bSJed Brown 
3439566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm,&localX));
344c4762a1bSJed Brown 
3459566063dSJacob Faibussowitsch   PetscCall(MatAssemblyBegin(J,MAT_FINAL_ASSEMBLY));
3469566063dSJacob Faibussowitsch   PetscCall(MatAssemblyEnd(J,MAT_FINAL_ASSEMBLY));
347c4762a1bSJed Brown   PetscFunctionReturn(0);
348c4762a1bSJed Brown }
349c4762a1bSJed Brown 
350c4762a1bSJed Brown PetscErrorCode SetInitialValues_Subnet(DM networkdm,Vec localX,PetscInt nv,PetscInt ne, const PetscInt *vtx, const PetscInt *edges,void* appctx)
351c4762a1bSJed Brown {
352c4762a1bSJed Brown   VERTEX_Power   bus;
353c4762a1bSJed Brown   PetscInt       i;
354c4762a1bSJed Brown   GEN            gen;
355c4762a1bSJed Brown   PetscBool      ghostvtex;
356c4762a1bSJed Brown   PetscScalar    *xarr;
357c4762a1bSJed Brown   PetscInt       key,numComps,j,offset;
358c4762a1bSJed Brown   void*          component;
359c4762a1bSJed Brown 
360c4762a1bSJed Brown   PetscFunctionBegin;
3619566063dSJacob Faibussowitsch   PetscCall(VecGetArray(localX,&xarr));
362c4762a1bSJed Brown   for (i = 0; i < nv; i++) {
3639566063dSJacob Faibussowitsch     PetscCall(DMNetworkIsGhostVertex(networkdm,vtx[i],&ghostvtex));
364c4762a1bSJed Brown     if (ghostvtex) continue;
365c4762a1bSJed Brown 
3669566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetLocalVecOffset(networkdm,vtx[i],ALL_COMPONENTS,&offset));
3679566063dSJacob Faibussowitsch     PetscCall(DMNetworkGetNumComponents(networkdm,vtx[i],&numComps));
368c4762a1bSJed Brown     for (j=0; j < numComps; j++) {
3699566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetComponent(networkdm,vtx[i],j,&key,&component,NULL));
370c4762a1bSJed Brown       if (key == 1) {
371c4762a1bSJed Brown         bus = (VERTEX_Power)(component);
372c4762a1bSJed Brown         xarr[offset] = bus->va*PETSC_PI/180.0;
373c4762a1bSJed Brown         xarr[offset+1] = bus->vm;
374c4762a1bSJed Brown       } else if (key == 2) {
375c4762a1bSJed Brown         gen = (GEN)(component);
376c4762a1bSJed Brown         if (!gen->status) continue;
377c4762a1bSJed Brown         xarr[offset+1] = gen->vs;
378c4762a1bSJed Brown         break;
379c4762a1bSJed Brown       }
380c4762a1bSJed Brown     }
381c4762a1bSJed Brown   }
3829566063dSJacob Faibussowitsch   PetscCall(VecRestoreArray(localX,&xarr));
383c4762a1bSJed Brown   PetscFunctionReturn(0);
384c4762a1bSJed Brown }
385c4762a1bSJed Brown 
386c4762a1bSJed Brown PetscErrorCode SetInitialValues(DM networkdm, Vec X,void* appctx)
387c4762a1bSJed Brown {
388c4762a1bSJed Brown   PetscInt       nv,ne;
389c4762a1bSJed Brown   const PetscInt *vtx,*edges;
390c4762a1bSJed Brown   Vec            localX;
391c4762a1bSJed Brown 
392c4762a1bSJed Brown   PetscFunctionBegin;
3939566063dSJacob Faibussowitsch   PetscCall(DMGetLocalVector(networkdm,&localX));
394c4762a1bSJed Brown 
3959566063dSJacob Faibussowitsch   PetscCall(VecSet(X,0.0));
3969566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalBegin(networkdm,X,INSERT_VALUES,localX));
3979566063dSJacob Faibussowitsch   PetscCall(DMGlobalToLocalEnd(networkdm,X,INSERT_VALUES,localX));
398c4762a1bSJed Brown 
399c4762a1bSJed Brown   /* Set initial guess for first subnetwork */
4009566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges));
4019566063dSJacob Faibussowitsch   PetscCall(SetInitialValues_Subnet(networkdm,localX,nv,ne,vtx,edges,appctx));
402c4762a1bSJed Brown 
403c4762a1bSJed Brown   /* Set initial guess for second subnetwork */
4049566063dSJacob Faibussowitsch   PetscCall(DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges));
4059566063dSJacob Faibussowitsch   PetscCall(SetInitialValues_Subnet(networkdm,localX,nv,ne,vtx,edges,appctx));
406c4762a1bSJed Brown 
4079566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalBegin(networkdm,localX,ADD_VALUES,X));
4089566063dSJacob Faibussowitsch   PetscCall(DMLocalToGlobalEnd(networkdm,localX,ADD_VALUES,X));
4099566063dSJacob Faibussowitsch   PetscCall(DMRestoreLocalVector(networkdm,&localX));
410c4762a1bSJed Brown   PetscFunctionReturn(0);
411c4762a1bSJed Brown }
412c4762a1bSJed Brown 
413c4762a1bSJed Brown int main(int argc,char ** argv)
414c4762a1bSJed Brown {
415c4762a1bSJed Brown   char             pfdata_file[PETSC_MAX_PATH_LEN]="case9.m";
416c4762a1bSJed Brown   PFDATA           *pfdata1,*pfdata2;
417f11a936eSBarry Smith   PetscInt         numEdges1=0,numEdges2=0;
4182bf73ac6SHong Zhang   PetscInt         *edgelist1 = NULL,*edgelist2 = NULL,componentkey[4];
419c4762a1bSJed Brown   DM               networkdm;
420c4762a1bSJed Brown   UserCtx_Power    User;
421956f8c0dSBarry Smith #if defined(PETSC_USE_LOG)
422c4762a1bSJed Brown   PetscLogStage    stage1,stage2;
423956f8c0dSBarry Smith #endif
424c4762a1bSJed Brown   PetscMPIInt      rank;
4252bf73ac6SHong Zhang   PetscInt         nsubnet = 2,nv,ne,i,j,genj,loadj;
4262bf73ac6SHong Zhang   const PetscInt   *vtx,*edges;
427c4762a1bSJed Brown   Vec              X,F;
428c4762a1bSJed Brown   Mat              J;
429c4762a1bSJed Brown   SNES             snes;
430c4762a1bSJed Brown 
431*327415f7SBarry Smith   PetscFunctionBeginUser;
4329566063dSJacob Faibussowitsch   PetscCall(PetscInitialize(&argc,&argv,"poweroptions",help));
4339566063dSJacob Faibussowitsch   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD,&rank));
434c4762a1bSJed Brown   {
435c4762a1bSJed 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 */
436c4762a1bSJed Brown     /* this is an experiment to see how the analyzer reacts */
437c4762a1bSJed Brown     const PetscMPIInt crank = rank;
438c4762a1bSJed Brown 
439c4762a1bSJed Brown     /* Create an empty network object */
4409566063dSJacob Faibussowitsch     PetscCall(DMNetworkCreate(PETSC_COMM_WORLD,&networkdm));
441c4762a1bSJed Brown 
442c4762a1bSJed Brown     /* Register the components in the network */
4439566063dSJacob Faibussowitsch     PetscCall(DMNetworkRegisterComponent(networkdm,"branchstruct",sizeof(struct _p_EDGE_Power),&componentkey[0]));
4449566063dSJacob Faibussowitsch     PetscCall(DMNetworkRegisterComponent(networkdm,"busstruct",sizeof(struct _p_VERTEX_Power),&componentkey[1]));
4459566063dSJacob Faibussowitsch     PetscCall(DMNetworkRegisterComponent(networkdm,"genstruct",sizeof(struct _p_GEN),&componentkey[2]));
4469566063dSJacob Faibussowitsch     PetscCall(DMNetworkRegisterComponent(networkdm,"loadstruct",sizeof(struct _p_LOAD),&componentkey[3]));
447c4762a1bSJed Brown 
4489566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister("Read Data",&stage1));
449c4762a1bSJed Brown     PetscLogStagePush(stage1);
450c4762a1bSJed Brown     /* READ THE DATA */
451c4762a1bSJed Brown     if (!crank) {
452c4762a1bSJed Brown       /* Only rank 0 reads the data */
4539566063dSJacob Faibussowitsch       PetscCall(PetscOptionsGetString(NULL,NULL,"-pfdata",pfdata_file,sizeof(pfdata_file),NULL));
454c4762a1bSJed Brown       /* HERE WE CREATE COPIES OF THE SAME NETWORK THAT WILL BE TREATED AS SUBNETWORKS */
455c4762a1bSJed Brown 
456c4762a1bSJed Brown       /*    READ DATA FOR THE FIRST SUBNETWORK */
4579566063dSJacob Faibussowitsch       PetscCall(PetscNew(&pfdata1));
4589566063dSJacob Faibussowitsch       PetscCall(PFReadMatPowerData(pfdata1,pfdata_file));
459c4762a1bSJed Brown       User.Sbase = pfdata1->sbase;
460c4762a1bSJed Brown 
461c4762a1bSJed Brown       numEdges1 = pfdata1->nbranch;
4629566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2*numEdges1,&edgelist1));
4639566063dSJacob Faibussowitsch       PetscCall(GetListofEdges_Power(pfdata1,edgelist1));
464c4762a1bSJed Brown 
465c4762a1bSJed Brown       /*    READ DATA FOR THE SECOND SUBNETWORK */
4669566063dSJacob Faibussowitsch       PetscCall(PetscNew(&pfdata2));
4679566063dSJacob Faibussowitsch       PetscCall(PFReadMatPowerData(pfdata2,pfdata_file));
468c4762a1bSJed Brown       User.Sbase = pfdata2->sbase;
469c4762a1bSJed Brown 
470c4762a1bSJed Brown       numEdges2 = pfdata2->nbranch;
4719566063dSJacob Faibussowitsch       PetscCall(PetscMalloc1(2*numEdges2,&edgelist2));
4729566063dSJacob Faibussowitsch       PetscCall(GetListofEdges_Power(pfdata2,edgelist2));
473c4762a1bSJed Brown     }
474c4762a1bSJed Brown 
475c4762a1bSJed Brown     PetscLogStagePop();
4769566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
4779566063dSJacob Faibussowitsch     PetscCall(PetscLogStageRegister("Create network",&stage2));
478c4762a1bSJed Brown     PetscLogStagePush(stage2);
479c4762a1bSJed Brown 
4802bf73ac6SHong Zhang     /* Set number of nodes/edges and edge connectivity */
4819566063dSJacob Faibussowitsch     PetscCall(DMNetworkSetNumSubNetworks(networkdm,PETSC_DECIDE,nsubnet));
4829566063dSJacob Faibussowitsch     PetscCall(DMNetworkAddSubnetwork(networkdm,"",numEdges1,edgelist1,NULL));
4839566063dSJacob Faibussowitsch     PetscCall(DMNetworkAddSubnetwork(networkdm,"",numEdges2,edgelist2,NULL));
484c4762a1bSJed Brown 
485c4762a1bSJed Brown     /* Set up the network layout */
4869566063dSJacob Faibussowitsch     PetscCall(DMNetworkLayoutSetUp(networkdm));
487c4762a1bSJed Brown 
488c4762a1bSJed Brown     /* Add network components only process 0 has any data to add*/
489c4762a1bSJed Brown     if (!crank) {
490c4762a1bSJed Brown       genj=0; loadj=0;
491c4762a1bSJed Brown 
492c4762a1bSJed Brown       /* ADD VARIABLES AND COMPONENTS FOR THE FIRST SUBNETWORK */
4939566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetSubnetwork(networkdm,0,&nv,&ne,&vtx,&edges));
494c4762a1bSJed Brown 
495c4762a1bSJed Brown       for (i = 0; i < ne; i++) {
4969566063dSJacob Faibussowitsch         PetscCall(DMNetworkAddComponent(networkdm,edges[i],componentkey[0],&pfdata1->branch[i],0));
497c4762a1bSJed Brown       }
498c4762a1bSJed Brown 
499c4762a1bSJed Brown       for (i = 0; i < nv; i++) {
5009566063dSJacob Faibussowitsch         PetscCall(DMNetworkAddComponent(networkdm,vtx[i],componentkey[1],&pfdata1->bus[i],2));
501c4762a1bSJed Brown         if (pfdata1->bus[i].ngen) {
502c4762a1bSJed Brown           for (j = 0; j < pfdata1->bus[i].ngen; j++) {
5039566063dSJacob Faibussowitsch             PetscCall(DMNetworkAddComponent(networkdm,vtx[i],componentkey[2],&pfdata1->gen[genj++],0));
504c4762a1bSJed Brown           }
505c4762a1bSJed Brown         }
506c4762a1bSJed Brown         if (pfdata1->bus[i].nload) {
507c4762a1bSJed Brown           for (j=0; j < pfdata1->bus[i].nload; j++) {
5089566063dSJacob Faibussowitsch             PetscCall(DMNetworkAddComponent(networkdm,vtx[i],componentkey[3],&pfdata1->load[loadj++],0));
509c4762a1bSJed Brown           }
510c4762a1bSJed Brown         }
511c4762a1bSJed Brown       }
512c4762a1bSJed Brown 
513c4762a1bSJed Brown       genj=0; loadj=0;
514c4762a1bSJed Brown 
515c4762a1bSJed Brown       /* ADD VARIABLES AND COMPONENTS FOR THE SECOND SUBNETWORK */
5169566063dSJacob Faibussowitsch       PetscCall(DMNetworkGetSubnetwork(networkdm,1,&nv,&ne,&vtx,&edges));
517c4762a1bSJed Brown 
518c4762a1bSJed Brown       for (i = 0; i < ne; i++) {
5199566063dSJacob Faibussowitsch         PetscCall(DMNetworkAddComponent(networkdm,edges[i],componentkey[0],&pfdata2->branch[i],0));
520c4762a1bSJed Brown       }
521c4762a1bSJed Brown 
522c4762a1bSJed Brown       for (i = 0; i < nv; i++) {
5239566063dSJacob Faibussowitsch         PetscCall(DMNetworkAddComponent(networkdm,vtx[i],componentkey[1],&pfdata2->bus[i],2));
524c4762a1bSJed Brown         if (pfdata2->bus[i].ngen) {
525c4762a1bSJed Brown           for (j = 0; j < pfdata2->bus[i].ngen; j++) {
5269566063dSJacob Faibussowitsch             PetscCall(DMNetworkAddComponent(networkdm,vtx[i],componentkey[2],&pfdata2->gen[genj++],0));
527c4762a1bSJed Brown           }
528c4762a1bSJed Brown         }
529c4762a1bSJed Brown         if (pfdata2->bus[i].nload) {
530c4762a1bSJed Brown           for (j=0; j < pfdata2->bus[i].nload; j++) {
5319566063dSJacob Faibussowitsch             PetscCall(DMNetworkAddComponent(networkdm,vtx[i],componentkey[3],&pfdata2->load[loadj++],0));
532c4762a1bSJed Brown           }
533c4762a1bSJed Brown         }
534c4762a1bSJed Brown       }
535c4762a1bSJed Brown     }
536c4762a1bSJed Brown 
537c4762a1bSJed Brown     /* Set up DM for use */
5389566063dSJacob Faibussowitsch     PetscCall(DMSetUp(networkdm));
539c4762a1bSJed Brown 
540c4762a1bSJed Brown     if (!crank) {
5419566063dSJacob Faibussowitsch       PetscCall(PetscFree(edgelist1));
5429566063dSJacob Faibussowitsch       PetscCall(PetscFree(edgelist2));
543c4762a1bSJed Brown     }
544c4762a1bSJed Brown 
545c4762a1bSJed Brown     if (!crank) {
5469566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata1->bus));
5479566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata1->gen));
5489566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata1->branch));
5499566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata1->load));
5509566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata1));
551c4762a1bSJed Brown 
5529566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata2->bus));
5539566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata2->gen));
5549566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata2->branch));
5559566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata2->load));
5569566063dSJacob Faibussowitsch       PetscCall(PetscFree(pfdata2));
557c4762a1bSJed Brown     }
558c4762a1bSJed Brown 
559c4762a1bSJed Brown     /* Distribute networkdm to multiple processes */
5609566063dSJacob Faibussowitsch     PetscCall(DMNetworkDistribute(&networkdm,0));
561c4762a1bSJed Brown 
562c4762a1bSJed Brown     PetscLogStagePop();
563c4762a1bSJed Brown 
564c4762a1bSJed Brown     /* Broadcast Sbase to all processors */
5659566063dSJacob Faibussowitsch     PetscCallMPI(MPI_Bcast(&User.Sbase,1,MPIU_SCALAR,0,PETSC_COMM_WORLD));
566c4762a1bSJed Brown 
5679566063dSJacob Faibussowitsch     PetscCall(DMCreateGlobalVector(networkdm,&X));
5689566063dSJacob Faibussowitsch     PetscCall(VecDuplicate(X,&F));
569c4762a1bSJed Brown 
5709566063dSJacob Faibussowitsch     PetscCall(DMCreateMatrix(networkdm,&J));
5719566063dSJacob Faibussowitsch     PetscCall(MatSetOption(J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_FALSE));
572c4762a1bSJed Brown 
5739566063dSJacob Faibussowitsch     PetscCall(SetInitialValues(networkdm,X,&User));
574c4762a1bSJed Brown 
575c4762a1bSJed Brown     /* HOOK UP SOLVER */
5769566063dSJacob Faibussowitsch     PetscCall(SNESCreate(PETSC_COMM_WORLD,&snes));
5779566063dSJacob Faibussowitsch     PetscCall(SNESSetDM(snes,networkdm));
5789566063dSJacob Faibussowitsch     PetscCall(SNESSetFunction(snes,F,FormFunction,&User));
5799566063dSJacob Faibussowitsch     PetscCall(SNESSetJacobian(snes,J,J,FormJacobian,&User));
5809566063dSJacob Faibussowitsch     PetscCall(SNESSetFromOptions(snes));
581c4762a1bSJed Brown 
5829566063dSJacob Faibussowitsch     PetscCall(SNESSolve(snes,NULL,X));
5839566063dSJacob Faibussowitsch     /* PetscCall(VecView(X,PETSC_VIEWER_STDOUT_WORLD)); */
584c4762a1bSJed Brown 
5859566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&X));
5869566063dSJacob Faibussowitsch     PetscCall(VecDestroy(&F));
5879566063dSJacob Faibussowitsch     PetscCall(MatDestroy(&J));
588c4762a1bSJed Brown 
5899566063dSJacob Faibussowitsch     PetscCall(SNESDestroy(&snes));
5909566063dSJacob Faibussowitsch     PetscCall(DMDestroy(&networkdm));
591c4762a1bSJed Brown   }
5929566063dSJacob Faibussowitsch   PetscCall(PetscFinalize());
593b122ec5aSJacob Faibussowitsch   return 0;
594c4762a1bSJed Brown }
595c4762a1bSJed Brown 
596c4762a1bSJed Brown /*TEST
597c4762a1bSJed Brown 
598c4762a1bSJed Brown    build:
599c4762a1bSJed Brown      depends: PFReadData.c pffunctions.c
600dfd57a17SPierre Jolivet      requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED)
601c4762a1bSJed Brown 
602c4762a1bSJed Brown    test:
603c4762a1bSJed Brown      args: -snes_rtol 1.e-3
604c4762a1bSJed Brown      localrunfiles: poweroptions case9.m
6052bf73ac6SHong Zhang      output_file: output/power_1.out
606c4762a1bSJed Brown 
607c4762a1bSJed Brown    test:
608c4762a1bSJed Brown      suffix: 2
609c4762a1bSJed Brown      args: -snes_rtol 1.e-3 -petscpartitioner_type simple
610c4762a1bSJed Brown      nsize: 4
611c4762a1bSJed Brown      localrunfiles: poweroptions case9.m
6122bf73ac6SHong Zhang      output_file: output/power_1.out
613c4762a1bSJed Brown 
614c4762a1bSJed Brown TEST*/
615