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