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