xref: /petsc/src/snes/tutorials/network/ex1.c (revision d16ceb7590f23315f034853cec344faff90b289e)
1 static char help[] = "This example demonstrates the use of DMNetwork with subnetworks for solving a coupled nonlinear \n\
2                       electric power grid and water pipe problem.\n\
3                       The available solver options are in the ex1options file \n\
4                       and the data files are in the datafiles of subdirectories.\n\
5                       Run this program: mpiexec -n <n> ./ex1 \n\\n";
6 
7 /*
8   Example:
9     mpiexec -n 3 ./ex1 -petscpartitioner_type parmetis -dmnetwork_view draw -dmnetwork_view_distributed draw -dmnetwork_view_rank_range 0,1,2
10 */
11 
12 #include "power/power.h"
13 #include "water/water.h"
14 
15 typedef struct {
16   UserCtx_Power appctx_power;
17   AppCtx_Water  appctx_water;
18   PetscInt      subsnes_id; /* snes solver id */
19   PetscInt      it;         /* iteration number */
20   Vec           localXold;  /* store previous solution, used by FormFunction_Dummy() */
21 } UserCtx;
22 
23 PetscErrorCode UserMonitor(SNES snes, PetscInt its, PetscReal fnorm, void *appctx)
24 {
25   UserCtx    *user = (UserCtx *)appctx;
26   Vec         X, localXold = user->localXold;
27   DM          networkdm;
28   PetscMPIInt rank;
29   MPI_Comm    comm;
30 
31   PetscFunctionBegin;
32   PetscCall(PetscObjectGetComm((PetscObject)snes, &comm));
33   PetscCallMPI(MPI_Comm_rank(comm, &rank));
34 #if 0
35   if (rank == 0) {
36     PetscInt       subsnes_id = user->subsnes_id;
37     if (subsnes_id == 2) {
38       PetscCall(PetscPrintf(PETSC_COMM_SELF," it %" PetscInt_FMT ", subsnes_id %" PetscInt_FMT ", fnorm %g\n",user->it,user->subsnes_id,(double)fnorm));
39     } else {
40       PetscCall(PetscPrintf(PETSC_COMM_SELF,"       subsnes_id %" PetscInt_FMT ", fnorm %g\n",user->subsnes_id,(double)fnorm));
41     }
42   }
43 #endif
44   PetscCall(SNESGetSolution(snes, &X));
45   PetscCall(SNESGetDM(snes, &networkdm));
46   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localXold));
47   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localXold));
48   PetscFunctionReturn(PETSC_SUCCESS);
49 }
50 
51 PetscErrorCode FormJacobian_subPower(SNES snes, Vec X, Mat J, Mat Jpre, void *appctx)
52 {
53   DM              networkdm;
54   Vec             localX;
55   PetscInt        nv, ne, i, j, offset, nvar, row;
56   const PetscInt *vtx, *edges;
57   PetscBool       ghostvtex;
58   PetscScalar     one = 1.0;
59   PetscMPIInt     rank;
60   MPI_Comm        comm;
61 
62   PetscFunctionBegin;
63   PetscCall(SNESGetDM(snes, &networkdm));
64   PetscCall(DMGetLocalVector(networkdm, &localX));
65 
66   PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm));
67   PetscCallMPI(MPI_Comm_rank(comm, &rank));
68 
69   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
70   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
71 
72   PetscCall(MatZeroEntries(J));
73 
74   /* Power subnetwork: copied from power/FormJacobian_Power() */
75   PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
76   PetscCall(FormJacobian_Power_private(networkdm, localX, J, nv, ne, vtx, edges, appctx));
77 
78   /* Water subnetwork: Identity */
79   PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges));
80   for (i = 0; i < nv; i++) {
81     PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex));
82     if (ghostvtex) continue;
83 
84     PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset));
85     PetscCall(DMNetworkGetComponent(networkdm, vtx[i], ALL_COMPONENTS, NULL, NULL, &nvar));
86     for (j = 0; j < nvar; j++) {
87       row = offset + j;
88       PetscCall(MatSetValues(J, 1, &row, 1, &row, &one, ADD_VALUES));
89     }
90   }
91   PetscCall(MatAssemblyBegin(J, MAT_FINAL_ASSEMBLY));
92   PetscCall(MatAssemblyEnd(J, MAT_FINAL_ASSEMBLY));
93 
94   PetscCall(DMRestoreLocalVector(networkdm, &localX));
95   PetscFunctionReturn(PETSC_SUCCESS);
96 }
97 
98 /* Dummy equation localF(X) = localX - localXold */
99 PetscErrorCode FormFunction_Dummy(DM networkdm, Vec localX, Vec localF, PetscInt nv, PetscInt ne, const PetscInt *vtx, const PetscInt *edges, void *appctx)
100 {
101   const PetscScalar *xarr, *xoldarr;
102   PetscScalar       *farr;
103   PetscInt           i, j, offset, nvar;
104   PetscBool          ghostvtex;
105   UserCtx           *user      = (UserCtx *)appctx;
106   Vec                localXold = user->localXold;
107 
108   PetscFunctionBegin;
109   PetscCall(VecGetArrayRead(localX, &xarr));
110   PetscCall(VecGetArrayRead(localXold, &xoldarr));
111   PetscCall(VecGetArray(localF, &farr));
112 
113   for (i = 0; i < nv; i++) {
114     PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghostvtex));
115     if (ghostvtex) continue;
116 
117     PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], ALL_COMPONENTS, &offset));
118     PetscCall(DMNetworkGetComponent(networkdm, vtx[i], ALL_COMPONENTS, NULL, NULL, &nvar));
119     for (j = 0; j < nvar; j++) farr[offset + j] = xarr[offset + j] - xoldarr[offset + j];
120   }
121 
122   PetscCall(VecRestoreArrayRead(localX, &xarr));
123   PetscCall(VecRestoreArrayRead(localXold, &xoldarr));
124   PetscCall(VecRestoreArray(localF, &farr));
125   PetscFunctionReturn(PETSC_SUCCESS);
126 }
127 
128 PetscErrorCode FormFunction(SNES snes, Vec X, Vec F, void *appctx)
129 {
130   DM              networkdm;
131   Vec             localX, localF;
132   PetscInt        nv, ne, v;
133   const PetscInt *vtx, *edges;
134   PetscMPIInt     rank;
135   MPI_Comm        comm;
136   UserCtx        *user         = (UserCtx *)appctx;
137   UserCtx_Power   appctx_power = (*user).appctx_power;
138   AppCtx_Water    appctx_water = (*user).appctx_water;
139 
140   PetscFunctionBegin;
141   PetscCall(SNESGetDM(snes, &networkdm));
142   PetscCall(PetscObjectGetComm((PetscObject)networkdm, &comm));
143   PetscCallMPI(MPI_Comm_rank(comm, &rank));
144 
145   PetscCall(DMGetLocalVector(networkdm, &localX));
146   PetscCall(DMGetLocalVector(networkdm, &localF));
147   PetscCall(VecSet(F, 0.0));
148   PetscCall(VecSet(localF, 0.0));
149 
150   PetscCall(DMGlobalToLocalBegin(networkdm, X, INSERT_VALUES, localX));
151   PetscCall(DMGlobalToLocalEnd(networkdm, X, INSERT_VALUES, localX));
152 
153   /* Form Function for power subnetwork */
154   PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
155   if (user->subsnes_id == 1) { /* snes_water only */
156     PetscCall(FormFunction_Dummy(networkdm, localX, localF, nv, ne, vtx, edges, user));
157   } else {
158     PetscCall(FormFunction_Power(networkdm, localX, localF, nv, ne, vtx, edges, &appctx_power));
159   }
160 
161   /* Form Function for water subnetwork */
162   PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges));
163   if (user->subsnes_id == 0) { /* snes_power only */
164     PetscCall(FormFunction_Dummy(networkdm, localX, localF, nv, ne, vtx, edges, user));
165   } else {
166     PetscCall(FormFunction_Water(networkdm, localX, localF, nv, ne, vtx, edges, NULL));
167   }
168 
169   /* Illustrate how to access the coupling vertex of the subnetworks without doing anything to F yet */
170   PetscCall(DMNetworkGetSharedVertices(networkdm, &nv, &vtx));
171   for (v = 0; v < nv; v++) {
172     PetscInt        key, ncomp, nvar, nconnedges, k, e, keye, goffset[3];
173     void           *component;
174     const PetscInt *connedges;
175 
176     PetscCall(DMNetworkGetComponent(networkdm, vtx[v], ALL_COMPONENTS, NULL, NULL, &nvar));
177     PetscCall(DMNetworkGetNumComponents(networkdm, vtx[v], &ncomp));
178     /* printf("  [%d] coupling vertex[%" PetscInt_FMT "]: v %" PetscInt_FMT ", ncomp %" PetscInt_FMT "; nvar %" PetscInt_FMT "\n",rank,v,vtx[v], ncomp,nvar); */
179 
180     for (k = 0; k < ncomp; k++) {
181       PetscCall(DMNetworkGetComponent(networkdm, vtx[v], k, &key, &component, &nvar));
182       PetscCall(DMNetworkGetGlobalVecOffset(networkdm, vtx[v], k, &goffset[k]));
183 
184       /* Verify the coupling vertex is a powernet load vertex or a water vertex */
185       switch (k) {
186       case 0:
187         PetscCheck(key == appctx_power.compkey_bus && nvar == 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "key %" PetscInt_FMT " not a power bus vertex or nvar %" PetscInt_FMT " != 2", key, nvar);
188         break;
189       case 1:
190         PetscCheck(key == appctx_power.compkey_load && nvar == 0 && goffset[1] == goffset[0] + 2, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not a power load vertex");
191         break;
192       case 2:
193         PetscCheck(key == appctx_water.compkey_vtx && nvar == 1 && goffset[2] == goffset[1], PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not a water vertex");
194         break;
195       default:
196         SETERRQ(PETSC_COMM_SELF, PETSC_ERR_ARG_OUTOFRANGE, "k %" PetscInt_FMT " is wrong", k);
197       }
198       /* printf("  [%d] coupling vertex[%" PetscInt_FMT "]: key %" PetscInt_FMT "; nvar %" PetscInt_FMT ", goffset %" PetscInt_FMT "\n",rank,v,key,nvar,goffset[k]); */
199     }
200 
201     /* Get its supporting edges */
202     PetscCall(DMNetworkGetSupportingEdges(networkdm, vtx[v], &nconnedges, &connedges));
203     /* printf("\n[%d] coupling vertex: nconnedges %" PetscInt_FMT "\n",rank,nconnedges); */
204     for (k = 0; k < nconnedges; k++) {
205       e = connedges[k];
206       PetscCall(DMNetworkGetNumComponents(networkdm, e, &ncomp));
207       /* printf("\n  [%d] connected edge[%" PetscInt_FMT "]=%" PetscInt_FMT " has ncomp %" PetscInt_FMT "\n",rank,k,e,ncomp); */
208       PetscCall(DMNetworkGetComponent(networkdm, e, 0, &keye, &component, NULL));
209       if (keye != appctx_water.compkey_edge) PetscCheck(keye == appctx_power.compkey_branch, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "Not a power branch");
210     }
211   }
212 
213   PetscCall(DMRestoreLocalVector(networkdm, &localX));
214 
215   PetscCall(DMLocalToGlobalBegin(networkdm, localF, ADD_VALUES, F));
216   PetscCall(DMLocalToGlobalEnd(networkdm, localF, ADD_VALUES, F));
217   PetscCall(DMRestoreLocalVector(networkdm, &localF));
218 #if 0
219   if (rank == 0) printf("F:\n");
220   PetscCall(VecView(F,PETSC_VIEWER_STDOUT_WORLD));
221 #endif
222   PetscFunctionReturn(PETSC_SUCCESS);
223 }
224 
225 PetscErrorCode SetInitialGuess(DM networkdm, Vec X, void *appctx)
226 {
227   PetscInt        nv, ne, i, j, ncomp, offset, key;
228   const PetscInt *vtx, *edges;
229   UserCtx        *user         = (UserCtx *)appctx;
230   Vec             localX       = user->localXold;
231   UserCtx_Power   appctx_power = (*user).appctx_power;
232   AppCtx_Water    appctx_water = (*user).appctx_water;
233   PetscBool       ghost;
234   PetscScalar    *xarr;
235   VERTEX_Power    bus;
236   VERTEX_Water    vertex;
237   void           *component;
238   GEN             gen;
239 
240   PetscFunctionBegin;
241   PetscCall(VecSet(X, 0.0));
242   PetscCall(VecSet(localX, 0.0));
243 
244   /* Set initial guess for power subnetwork */
245   PetscCall(DMNetworkGetSubnetwork(networkdm, 0, &nv, &ne, &vtx, &edges));
246   PetscCall(SetInitialGuess_Power(networkdm, localX, nv, ne, vtx, edges, &appctx_power));
247 
248   /* Set initial guess for water subnetwork */
249   PetscCall(DMNetworkGetSubnetwork(networkdm, 1, &nv, &ne, &vtx, &edges));
250   PetscCall(SetInitialGuess_Water(networkdm, localX, nv, ne, vtx, edges, NULL));
251 
252   /* Set initial guess at the coupling vertex */
253   PetscCall(VecGetArray(localX, &xarr));
254   PetscCall(DMNetworkGetSharedVertices(networkdm, &nv, &vtx));
255   for (i = 0; i < nv; i++) {
256     PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[i], &ghost));
257     if (ghost) continue;
258 
259     PetscCall(DMNetworkGetNumComponents(networkdm, vtx[i], &ncomp));
260     for (j = 0; j < ncomp; j++) {
261       PetscCall(DMNetworkGetLocalVecOffset(networkdm, vtx[i], j, &offset));
262       PetscCall(DMNetworkGetComponent(networkdm, vtx[i], j, &key, (void **)&component, NULL));
263       if (key == appctx_power.compkey_bus) {
264         bus              = (VERTEX_Power)(component);
265         xarr[offset]     = bus->va * PETSC_PI / 180.0;
266         xarr[offset + 1] = bus->vm;
267       } else if (key == appctx_power.compkey_gen) {
268         gen = (GEN)(component);
269         if (!gen->status) continue;
270         xarr[offset + 1] = gen->vs;
271       } else if (key == appctx_water.compkey_vtx) {
272         vertex = (VERTEX_Water)(component);
273         if (vertex->type == VERTEX_TYPE_JUNCTION) {
274           xarr[offset] = 100;
275         } else if (vertex->type == VERTEX_TYPE_RESERVOIR) {
276           xarr[offset] = vertex->res.head;
277         } else {
278           xarr[offset] = vertex->tank.initlvl + vertex->tank.elev;
279         }
280       }
281     }
282   }
283   PetscCall(VecRestoreArray(localX, &xarr));
284 
285   PetscCall(DMLocalToGlobalBegin(networkdm, localX, ADD_VALUES, X));
286   PetscCall(DMLocalToGlobalEnd(networkdm, localX, ADD_VALUES, X));
287   PetscFunctionReturn(PETSC_SUCCESS);
288 }
289 
290 /* Set coordinates */
291 static PetscErrorCode CoordinateVecSetUp(DM dm, Vec coords)
292 {
293   DM              dmclone;
294   PetscInt        i, gidx, offset, v, nv, Nsubnet;
295   const PetscInt *vtx;
296   PetscScalar    *carray;
297 
298   PetscFunctionBeginUser;
299   PetscCall(DMGetCoordinateDM(dm, &dmclone));
300   PetscCall(VecGetArrayWrite(coords, &carray));
301   PetscCall(DMNetworkGetNumSubNetworks(dm, NULL, &Nsubnet));
302   for (i = 0; i < Nsubnet; i++) {
303     PetscCall(DMNetworkGetSubnetwork(dm, i, &nv, NULL, &vtx, NULL));
304     for (v = 0; v < nv; v++) {
305       PetscCall(DMNetworkGetGlobalVertexIndex(dm, vtx[v], &gidx));
306       PetscCall(DMNetworkGetLocalVecOffset(dmclone, vtx[v], 0, &offset));
307       switch (gidx) {
308       case 0:
309         carray[offset]     = -1.0;
310         carray[offset + 1] = -1.0;
311         break;
312       case 1:
313         carray[offset]     = -2.0;
314         carray[offset + 1] = 2.0;
315         break;
316       case 2:
317         carray[offset]     = 0.0;
318         carray[offset + 1] = 2.0;
319         break;
320       case 3:
321         carray[offset]     = -1.0;
322         carray[offset + 1] = 0.0;
323         break;
324       case 4:
325         carray[offset]     = 0.0;
326         carray[offset + 1] = 0.0;
327         break;
328       case 5:
329         carray[offset]     = 0.0;
330         carray[offset + 1] = 1.0;
331         break;
332       case 6:
333         carray[offset]     = -1.0;
334         carray[offset + 1] = 1.0;
335         break;
336       case 7:
337         carray[offset]     = -2.0;
338         carray[offset + 1] = 1.0;
339         break;
340       case 8:
341         carray[offset]     = -2.0;
342         carray[offset + 1] = 0.0;
343         break;
344       case 9:
345         carray[offset]     = 1.0;
346         carray[offset + 1] = 0.0;
347         break;
348       case 10:
349         carray[offset]     = 1.0;
350         carray[offset + 1] = -1.0;
351         break;
352       case 11:
353         carray[offset]     = 2.0;
354         carray[offset + 1] = -1.0;
355         break;
356       case 12:
357         carray[offset]     = 2.0;
358         carray[offset + 1] = 0.0;
359         break;
360       case 13:
361         carray[offset]     = 0.0;
362         carray[offset + 1] = -1.0;
363         break;
364       case 14:
365         carray[offset]     = 2.0;
366         carray[offset + 1] = 1.0;
367         break;
368       default:
369         PetscCheck(gidx < 15 && gidx > -1, PETSC_COMM_SELF, PETSC_ERR_ARG_WRONG, "gidx %" PetscInt_FMT "must between 0 and 14", gidx);
370       }
371     }
372   }
373   PetscCall(VecRestoreArrayWrite(coords, &carray));
374   PetscFunctionReturn(PETSC_SUCCESS);
375 }
376 
377 static PetscErrorCode CoordinatePrint(DM dm)
378 {
379   DM                 dmclone;
380   PetscInt           cdim, v, off, vglobal, vStart, vEnd;
381   const PetscScalar *carray;
382   Vec                coords;
383   MPI_Comm           comm;
384   PetscMPIInt        rank;
385 
386   PetscFunctionBegin;
387   PetscCall(PetscObjectGetComm((PetscObject)dm, &comm));
388   PetscCallMPI(MPI_Comm_rank(comm, &rank));
389 
390   PetscCall(DMGetCoordinateDM(dm, &dmclone));
391   PetscCall(DMNetworkGetVertexRange(dm, &vStart, &vEnd));
392   PetscCall(DMGetCoordinatesLocal(dm, &coords));
393 
394   PetscCall(DMGetCoordinateDim(dm, &cdim));
395   PetscCall(VecGetArrayRead(coords, &carray));
396 
397   PetscCall(PetscPrintf(MPI_COMM_WORLD, "\nCoordinatePrint, cdim %" PetscInt_FMT ":\n", cdim));
398   PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "[%i]\n", rank));
399   for (v = vStart; v < vEnd; v++) {
400     PetscCall(DMNetworkGetLocalVecOffset(dmclone, v, 0, &off));
401     PetscCall(DMNetworkGetGlobalVertexIndex(dmclone, v, &vglobal));
402     switch (cdim) {
403     case 2:
404       PetscCall(PetscSynchronizedPrintf(MPI_COMM_WORLD, "Vertex: %" PetscInt_FMT ", x =  %f y = %f \n", vglobal, (double)PetscRealPart(carray[off]), (double)PetscRealPart(carray[off + 1])));
405       break;
406     default:
407       PetscCheck(cdim == 2, MPI_COMM_WORLD, PETSC_ERR_SUP, "Only supports Network embedding dimension of 2, not supplied  %" PetscInt_FMT, cdim);
408       break;
409     }
410   }
411   PetscCall(PetscSynchronizedFlush(MPI_COMM_WORLD, NULL));
412   PetscCall(VecRestoreArrayRead(coords, &carray));
413   PetscFunctionReturn(PETSC_SUCCESS);
414 }
415 
416 int main(int argc, char **argv)
417 {
418   DM                  networkdm, dmclone;
419   PetscMPIInt         rank, size;
420   PetscInt            Nsubnet = 2, numVertices[2], numEdges[2], i, j, nv, ne, it_max = 10;
421   PetscInt            vStart, vEnd, compkey;
422   const PetscInt     *vtx, *edges;
423   Vec                 X, F, coords;
424   SNES                snes, snes_power, snes_water;
425   Mat                 Jac;
426   PetscBool           ghost, viewJ = PETSC_FALSE, viewX = PETSC_FALSE, test = PETSC_FALSE, distribute = PETSC_TRUE, flg, printCoord = PETSC_FALSE, viewCSV = PETSC_FALSE;
427   UserCtx             user;
428   SNESConvergedReason reason;
429 #if defined(PETSC_USE_LOG)
430   PetscLogStage stage[4];
431 #endif
432 
433   /* Power subnetwork */
434   UserCtx_Power *appctx_power                    = &user.appctx_power;
435   char           pfdata_file[PETSC_MAX_PATH_LEN] = "power/case9.m";
436   PFDATA        *pfdata                          = NULL;
437   PetscInt       genj, loadj, *edgelist_power = NULL, power_netnum;
438   PetscScalar    Sbase = 0.0;
439 
440   /* Water subnetwork */
441   AppCtx_Water *appctx_water                       = &user.appctx_water;
442   WATERDATA    *waterdata                          = NULL;
443   char          waterdata_file[PETSC_MAX_PATH_LEN] = "water/sample1.inp";
444   PetscInt     *edgelist_water                     = NULL, water_netnum;
445 
446   /* Shared vertices between subnetworks */
447   PetscInt power_svtx, water_svtx;
448 
449   PetscFunctionBeginUser;
450   PetscCall(PetscInitialize(&argc, &argv, "ex1options", help));
451   PetscCallMPI(MPI_Comm_rank(PETSC_COMM_WORLD, &rank));
452   PetscCallMPI(MPI_Comm_size(PETSC_COMM_WORLD, &size));
453 
454   /* (1) Read Data - Only rank 0 reads the data */
455   PetscCall(PetscLogStageRegister("Read Data", &stage[0]));
456   PetscCall(PetscLogStagePush(stage[0]));
457 
458   for (i = 0; i < Nsubnet; i++) {
459     numVertices[i] = 0;
460     numEdges[i]    = 0;
461   }
462 
463   /* All processes READ THE DATA FOR THE FIRST SUBNETWORK: Electric Power Grid */
464   /* Used for shared vertex, because currently the coupling info must be available in all processes!!! */
465   PetscCall(PetscOptionsGetString(NULL, NULL, "-pfdata", pfdata_file, PETSC_MAX_PATH_LEN - 1, NULL));
466   PetscCall(PetscNew(&pfdata));
467   PetscCall(PFReadMatPowerData(pfdata, pfdata_file));
468   if (rank == 0) PetscCall(PetscPrintf(PETSC_COMM_SELF, "Power network: nb = %" PetscInt_FMT ", ngen = %" PetscInt_FMT ", nload = %" PetscInt_FMT ", nbranch = %" PetscInt_FMT "\n", pfdata->nbus, pfdata->ngen, pfdata->nload, pfdata->nbranch));
469   Sbase = pfdata->sbase;
470   if (rank == 0) { /* proc[0] will create Electric Power Grid */
471     numEdges[0]    = pfdata->nbranch;
472     numVertices[0] = pfdata->nbus;
473 
474     PetscCall(PetscMalloc1(2 * numEdges[0], &edgelist_power));
475     PetscCall(GetListofEdges_Power(pfdata, edgelist_power));
476   }
477   /* Broadcast power Sbase to all processors */
478   PetscCallMPI(MPI_Bcast(&Sbase, 1, MPIU_SCALAR, 0, PETSC_COMM_WORLD));
479   appctx_power->Sbase     = Sbase;
480   appctx_power->jac_error = PETSC_FALSE;
481   /* If external option activated. Introduce error in jacobian */
482   PetscCall(PetscOptionsHasName(NULL, NULL, "-jac_error", &appctx_power->jac_error));
483 
484   /* All processes READ THE DATA FOR THE SECOND SUBNETWORK: Water */
485   /* Used for shared vertex, because currently the coupling info must be available in all processes!!! */
486   PetscCall(PetscNew(&waterdata));
487   PetscCall(PetscOptionsGetString(NULL, NULL, "-waterdata", waterdata_file, PETSC_MAX_PATH_LEN - 1, NULL));
488   PetscCall(WaterReadData(waterdata, waterdata_file));
489   if (size == 1 || (size > 1 && rank == 1)) {
490     PetscCall(PetscCalloc1(2 * waterdata->nedge, &edgelist_water));
491     PetscCall(GetListofEdges_Water(waterdata, edgelist_water));
492     numEdges[1]    = waterdata->nedge;
493     numVertices[1] = waterdata->nvertex;
494   }
495   PetscCall(PetscLogStagePop());
496 
497   /* (2) Create a network consist of two subnetworks */
498   PetscCall(PetscLogStageRegister("Net Setup", &stage[1]));
499   PetscCall(PetscLogStagePush(stage[1]));
500 
501   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewCSV", &viewCSV, NULL));
502   PetscCall(PetscOptionsGetBool(NULL, NULL, "-printCoord", &printCoord, NULL));
503   PetscCall(PetscOptionsGetBool(NULL, NULL, "-test", &test, NULL));
504   PetscCall(PetscOptionsGetBool(NULL, NULL, "-distribute", &distribute, NULL));
505 
506   /* Create an empty network object */
507   PetscCall(DMNetworkCreate(PETSC_COMM_WORLD, &networkdm));
508 
509   /* Register the components in the network */
510   PetscCall(DMNetworkRegisterComponent(networkdm, "branchstruct", sizeof(struct _p_EDGE_Power), &appctx_power->compkey_branch));
511   PetscCall(DMNetworkRegisterComponent(networkdm, "busstruct", sizeof(struct _p_VERTEX_Power), &appctx_power->compkey_bus));
512   PetscCall(DMNetworkRegisterComponent(networkdm, "genstruct", sizeof(struct _p_GEN), &appctx_power->compkey_gen));
513   PetscCall(DMNetworkRegisterComponent(networkdm, "loadstruct", sizeof(struct _p_LOAD), &appctx_power->compkey_load));
514 
515   PetscCall(DMNetworkRegisterComponent(networkdm, "edge_water", sizeof(struct _p_EDGE_Water), &appctx_water->compkey_edge));
516   PetscCall(DMNetworkRegisterComponent(networkdm, "vertex_water", sizeof(struct _p_VERTEX_Water), &appctx_water->compkey_vtx));
517 
518   PetscCall(PetscSynchronizedPrintf(PETSC_COMM_WORLD, "[%d] Total local nvertices %" PetscInt_FMT " + %" PetscInt_FMT " = %" PetscInt_FMT ", nedges %" PetscInt_FMT " + %" PetscInt_FMT " = %" PetscInt_FMT "\n", rank, numVertices[0], numVertices[1], numVertices[0] + numVertices[1], numEdges[0], numEdges[1], numEdges[0] + numEdges[1]));
519   PetscCall(PetscSynchronizedFlush(PETSC_COMM_WORLD, PETSC_STDOUT));
520 
521   PetscCall(DMNetworkSetNumSubNetworks(networkdm, PETSC_DECIDE, Nsubnet));
522   PetscCall(DMNetworkAddSubnetwork(networkdm, "power", numEdges[0], edgelist_power, &power_netnum));
523   PetscCall(DMNetworkAddSubnetwork(networkdm, "water", numEdges[1], edgelist_water, &water_netnum));
524 
525   /* vertex subnet[0].4 shares with vertex subnet[1].0 */
526   power_svtx = 4;
527   water_svtx = 0;
528   PetscCall(DMNetworkAddSharedVertices(networkdm, power_netnum, water_netnum, 1, &power_svtx, &water_svtx));
529 
530   /* Set up the network layout */
531   PetscCall(DMNetworkLayoutSetUp(networkdm));
532 
533   /* ADD VARIABLES AND COMPONENTS FOR THE POWER SUBNETWORK */
534   genj  = 0;
535   loadj = 0;
536   PetscCall(DMNetworkGetSubnetwork(networkdm, power_netnum, &nv, &ne, &vtx, &edges));
537 
538   for (i = 0; i < ne; i++) PetscCall(DMNetworkAddComponent(networkdm, edges[i], appctx_power->compkey_branch, &pfdata->branch[i], 0));
539 
540   for (i = 0; i < nv; i++) {
541     PetscCall(DMNetworkIsSharedVertex(networkdm, vtx[i], &flg));
542     if (flg) continue;
543 
544     PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_bus, &pfdata->bus[i], 2));
545     if (pfdata->bus[i].ngen) {
546       for (j = 0; j < pfdata->bus[i].ngen; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_gen, &pfdata->gen[genj++], 0));
547     }
548     if (pfdata->bus[i].nload) {
549       for (j = 0; j < pfdata->bus[i].nload; j++) PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_load, &pfdata->load[loadj++], 0));
550     }
551   }
552 
553   /* ADD VARIABLES AND COMPONENTS FOR THE WATER SUBNETWORK */
554   PetscCall(DMNetworkGetSubnetwork(networkdm, water_netnum, &nv, &ne, &vtx, &edges));
555   for (i = 0; i < ne; i++) PetscCall(DMNetworkAddComponent(networkdm, edges[i], appctx_water->compkey_edge, &waterdata->edge[i], 0));
556 
557   for (i = 0; i < nv; i++) {
558     PetscCall(DMNetworkIsSharedVertex(networkdm, vtx[i], &flg));
559     if (flg) continue;
560 
561     PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_water->compkey_vtx, &waterdata->vertex[i], 1));
562   }
563 
564   /* ADD VARIABLES AND COMPONENTS AT THE SHARED VERTEX: net[0].4 coupls with net[1].0 -- owning and all ghost ranks of the vertex do this */
565   PetscCall(DMNetworkGetSharedVertices(networkdm, &nv, &vtx));
566   for (i = 0; i < nv; i++) {
567     /* power */
568     PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_bus, &pfdata->bus[4], 2));
569     /* bus[4] is a load, add its component */
570     PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_power->compkey_load, &pfdata->load[0], 0));
571 
572     /* water */
573     PetscCall(DMNetworkAddComponent(networkdm, vtx[i], appctx_water->compkey_vtx, &waterdata->vertex[0], 1));
574   }
575 
576   /* Set coordinates for visualization */
577   PetscCall(DMSetCoordinateDim(networkdm, 2));
578   PetscCall(DMGetCoordinateDM(networkdm, &dmclone));
579   PetscCall(DMNetworkGetVertexRange(dmclone, &vStart, &vEnd));
580   PetscCall(DMNetworkRegisterComponent(dmclone, "coordinates", 0, &compkey));
581   for (i = vStart; i < vEnd; i++) PetscCall(DMNetworkAddComponent(dmclone, i, compkey, NULL, 2));
582   PetscCall(DMNetworkFinalizeComponents(dmclone));
583 
584   PetscCall(DMCreateLocalVector(dmclone, &coords));
585   PetscCall(DMSetCoordinatesLocal(networkdm, coords)); /* set/get coords to/from networkdm */
586   PetscCall(CoordinateVecSetUp(networkdm, coords));
587   if (printCoord) PetscCall(CoordinatePrint(networkdm));
588 
589   /* Set up DM for use */
590   PetscCall(DMSetFromOptions(networkdm));
591   PetscCall(DMSetUp(networkdm));
592 
593   /* Free user objects */
594   PetscCall(PetscFree(edgelist_power));
595   PetscCall(PetscFree(pfdata->bus));
596   PetscCall(PetscFree(pfdata->gen));
597   PetscCall(PetscFree(pfdata->branch));
598   PetscCall(PetscFree(pfdata->load));
599   PetscCall(PetscFree(pfdata));
600 
601   PetscCall(PetscFree(edgelist_water));
602   PetscCall(PetscFree(waterdata->vertex));
603   PetscCall(PetscFree(waterdata->edge));
604   PetscCall(PetscFree(waterdata));
605 
606   /* Re-distribute networkdm to multiple processes for better job balance */
607   if (distribute) {
608     PetscCall(DMNetworkDistribute(&networkdm, 0));
609 
610     if (printCoord) PetscCall(CoordinatePrint(networkdm));
611     if (viewCSV) { /* CSV View of network with coordinates */
612       PetscCall(PetscViewerPushFormat(PETSC_VIEWER_STDOUT_WORLD, PETSC_VIEWER_ASCII_CSV));
613       PetscCall(DMView(networkdm, PETSC_VIEWER_STDOUT_WORLD));
614       PetscCall(PetscViewerPopFormat(PETSC_VIEWER_STDOUT_WORLD));
615     }
616   }
617   PetscCall(VecDestroy(&coords));
618 
619   /* Test DMNetworkGetSubnetwork() and DMNetworkGetSubnetworkSharedVertices() */
620   if (test) {
621     PetscInt v, gidx;
622     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
623     for (i = 0; i < Nsubnet; i++) {
624       PetscCall(DMNetworkGetSubnetwork(networkdm, i, &nv, &ne, &vtx, &edges));
625       PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] After distribute, subnet[%" PetscInt_FMT "] ne %" PetscInt_FMT ", nv %" PetscInt_FMT "\n", rank, i, ne, nv));
626       PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
627 
628       for (v = 0; v < nv; v++) {
629         PetscCall(DMNetworkIsGhostVertex(networkdm, vtx[v], &ghost));
630         PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, vtx[v], &gidx));
631         PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] subnet[%" PetscInt_FMT "] v %" PetscInt_FMT " %" PetscInt_FMT "; ghost %d\n", rank, i, vtx[v], gidx, ghost));
632       }
633       PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
634     }
635     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
636 
637     PetscCall(DMNetworkGetSharedVertices(networkdm, &nv, &vtx));
638     PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] After distribute, num of shared vertices nsv = %" PetscInt_FMT "\n", rank, nv));
639     for (v = 0; v < nv; v++) {
640       PetscCall(DMNetworkGetGlobalVertexIndex(networkdm, vtx[v], &gidx));
641       PetscCall(PetscPrintf(PETSC_COMM_SELF, "[%d] sv %" PetscInt_FMT ", gidx=%" PetscInt_FMT "\n", rank, vtx[v], gidx));
642     }
643     PetscCallMPI(MPI_Barrier(PETSC_COMM_WORLD));
644   }
645 
646   /* Create solution vector X */
647   PetscCall(DMCreateGlobalVector(networkdm, &X));
648   PetscCall(VecDuplicate(X, &F));
649   PetscCall(DMGetLocalVector(networkdm, &user.localXold));
650   PetscCall(PetscLogStagePop());
651 
652   /* (3) Setup Solvers */
653   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewJ", &viewJ, NULL));
654   PetscCall(PetscOptionsGetBool(NULL, NULL, "-viewX", &viewX, NULL));
655 
656   PetscCall(PetscLogStageRegister("SNES Setup", &stage[2]));
657   PetscCall(PetscLogStagePush(stage[2]));
658 
659   PetscCall(SetInitialGuess(networkdm, X, &user));
660 
661   /* Create coupled snes */
662   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES_coupled setup ......\n"));
663   user.subsnes_id = Nsubnet;
664   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes));
665   PetscCall(SNESSetDM(snes, networkdm));
666   PetscCall(SNESSetOptionsPrefix(snes, "coupled_"));
667   PetscCall(SNESSetFunction(snes, F, FormFunction, &user));
668   PetscCall(SNESMonitorSet(snes, UserMonitor, &user, NULL));
669   PetscCall(SNESSetFromOptions(snes));
670 
671   if (viewJ) {
672     /* View Jac structure */
673     PetscCall(SNESGetJacobian(snes, &Jac, NULL, NULL, NULL));
674     PetscCall(MatView(Jac, PETSC_VIEWER_DRAW_WORLD));
675   }
676 
677   if (viewX) {
678     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Solution:\n"));
679     PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
680   }
681 
682   if (viewJ) {
683     /* View assembled Jac */
684     PetscCall(MatView(Jac, PETSC_VIEWER_DRAW_WORLD));
685   }
686 
687   /* Create snes_power */
688   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES_power setup ......\n"));
689 
690   user.subsnes_id = 0;
691   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_power));
692   PetscCall(SNESSetDM(snes_power, networkdm));
693   PetscCall(SNESSetOptionsPrefix(snes_power, "power_"));
694   PetscCall(SNESSetFunction(snes_power, F, FormFunction, &user));
695   PetscCall(SNESMonitorSet(snes_power, UserMonitor, &user, NULL));
696 
697   /* Use user-provide Jacobian */
698   PetscCall(DMCreateMatrix(networkdm, &Jac));
699   PetscCall(SNESSetJacobian(snes_power, Jac, Jac, FormJacobian_subPower, &user));
700   PetscCall(SNESSetFromOptions(snes_power));
701 
702   if (viewX) {
703     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Power Solution:\n"));
704     PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
705   }
706   if (viewJ) {
707     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Power Jac:\n"));
708     PetscCall(SNESGetJacobian(snes_power, &Jac, NULL, NULL, NULL));
709     PetscCall(MatView(Jac, PETSC_VIEWER_DRAW_WORLD));
710     /* PetscCall(MatView(Jac,PETSC_VIEWER_STDOUT_WORLD)); */
711   }
712 
713   /* Create snes_water */
714   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "SNES_water setup......\n"));
715 
716   user.subsnes_id = 1;
717   PetscCall(SNESCreate(PETSC_COMM_WORLD, &snes_water));
718   PetscCall(SNESSetDM(snes_water, networkdm));
719   PetscCall(SNESSetOptionsPrefix(snes_water, "water_"));
720   PetscCall(SNESSetFunction(snes_water, F, FormFunction, &user));
721   PetscCall(SNESMonitorSet(snes_water, UserMonitor, &user, NULL));
722   PetscCall(SNESSetFromOptions(snes_water));
723 
724   if (viewX) {
725     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Water Solution:\n"));
726     PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
727   }
728   PetscCall(PetscLogStagePop());
729 
730   /* (4) Solve */
731   PetscCall(PetscLogStageRegister("SNES Solve", &stage[3]));
732   PetscCall(PetscLogStagePush(stage[3]));
733   user.it = 0;
734   reason  = SNES_DIVERGED_DTOL;
735   while (user.it < it_max && (PetscInt)reason < 0) {
736 #if 0
737     user.subsnes_id = 0;
738     PetscCall(SNESSolve(snes_power,NULL,X));
739 
740     user.subsnes_id = 1;
741     PetscCall(SNESSolve(snes_water,NULL,X));
742 #endif
743     user.subsnes_id = Nsubnet;
744     PetscCall(SNESSolve(snes, NULL, X));
745 
746     PetscCall(SNESGetConvergedReason(snes, &reason));
747     user.it++;
748   }
749   PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Coupled_SNES converged in %" PetscInt_FMT " iterations\n", user.it));
750   if (viewX) {
751     PetscCall(PetscPrintf(PETSC_COMM_WORLD, "Final Solution:\n"));
752     PetscCall(VecView(X, PETSC_VIEWER_STDOUT_WORLD));
753   }
754   PetscCall(PetscLogStagePop());
755 
756   /* Free objects */
757   /* -------------*/
758   PetscCall(VecDestroy(&X));
759   PetscCall(VecDestroy(&F));
760   PetscCall(DMRestoreLocalVector(networkdm, &user.localXold));
761 
762   PetscCall(SNESDestroy(&snes));
763   PetscCall(MatDestroy(&Jac));
764   PetscCall(SNESDestroy(&snes_power));
765   PetscCall(SNESDestroy(&snes_water));
766 
767   PetscCall(DMDestroy(&networkdm));
768   PetscCall(PetscFinalize());
769   return 0;
770 }
771 
772 /*TEST
773 
774    build:
775      requires: !complex double defined(PETSC_HAVE_ATTRIBUTEALIGNED)
776      depends: power/PFReadData.c power/pffunctions.c water/waterreaddata.c water/waterfunctions.c
777 
778    test:
779       args: -coupled_snes_converged_reason -options_left no -dmnetwork_view -fp_trap 0
780       localrunfiles: ex1options power/case9.m water/sample1.inp
781       output_file: output/ex1.out
782 
783    test:
784       suffix: 2
785       nsize: 3
786       args: -coupled_snes_converged_reason -options_left no -petscpartitioner_type parmetis -fp_trap 0
787       localrunfiles: ex1options power/case9.m water/sample1.inp
788       output_file: output/ex1_2.out
789       requires: parmetis
790 
791    test:
792       suffix: 3
793       nsize: 3
794       args: -coupled_snes_converged_reason -options_left no -distribute false -fp_trap 0
795       localrunfiles: ex1options power/case9.m water/sample1.inp
796       output_file: output/ex1_2.out
797 
798    test:
799       suffix: 4
800       nsize: 4
801       args: -coupled_snes_converged_reason -options_left no -petscpartitioner_type simple -dmnetwork_view -dmnetwork_view_distributed -fp_trap 0
802       localrunfiles: ex1options power/case9.m water/sample1.inp
803       output_file: output/ex1_4.out
804 
805    test:
806       suffix: 5
807       args: -coupled_snes_converged_reason -options_left no -viewCSV -fp_trap 0
808       localrunfiles: ex1options power/case9.m water/sample1.inp
809       output_file: output/ex1_5.out
810 
811    test:
812       suffix: 6
813       nsize: 3
814       args: -coupled_snes_converged_reason -options_left no -petscpartitioner_type parmetis -dmnetwork_view_distributed draw:null -fp_trap 0
815       localrunfiles: ex1options power/case9.m water/sample1.inp
816       output_file: output/ex1_2.out
817       requires: parmetis
818 
819 TEST*/
820