xref: /petsc/src/dm/impls/network/network.c (revision a6f8f0eb0f15cce44058a6ad44e802f69a95fee6)
1 #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
2 
3 /*
4   Creates the component header and value objects for a network point
5 */
6 static PetscErrorCode SetUpNetworkHeaderComponentValue(DM dm,DMNetworkComponentHeader header,DMNetworkComponentValue cvalue)
7 {
8   PetscFunctionBegin;
9   /* Allocate arrays for component information */
10   PetscCall(PetscCalloc5(header->maxcomps,&header->size,header->maxcomps,&header->key,header->maxcomps,&header->offset,header->maxcomps,&header->nvar,header->maxcomps,&header->offsetvarrel));
11   PetscCall(PetscCalloc1(header->maxcomps,&cvalue->data));
12 
13   /* The size of the header is the size of struct _p_DMNetworkComponentHeader. Since the struct contains PetscInt pointers we cannot use sizeof(struct). So, we need to explicitly calculate the size.
14    If the data header struct changes then this header size calculation needs to be updated. */
15   header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt);
16   header->hsize /= sizeof(DMNetworkComponentGenericDataType);
17   PetscFunctionReturn(0);
18 }
19 
20 /*@
21   DMNetworkGetPlex - Gets the Plex DM associated with this network DM
22 
23   Not collective
24 
25   Input Parameters:
26 . dm - the dm object
27 
28   Output Parameters:
29 . plexdm - the plex dm object
30 
31   Level: Advanced
32 
33 .seealso: DMNetworkCreate()
34 @*/
35 PetscErrorCode DMNetworkGetPlex(DM dm,DM *plexdm)
36 {
37   DM_Network *network = (DM_Network*)dm->data;
38 
39   PetscFunctionBegin;
40   *plexdm = network->plex;
41   PetscFunctionReturn(0);
42 }
43 
44 /*@
45   DMNetworkGetNumSubNetworks - Gets the the number of subnetworks
46 
47   Not collective
48 
49   Input Parameter:
50 . dm - the dm object
51 
52   Output Parameters:
53 + nsubnet - local number of subnetworks
54 - Nsubnet - global number of subnetworks
55 
56   Level: beginner
57 
58 .seealso: DMNetworkCreate(), DMNetworkSetNumSubNetworks()
59 @*/
60 PetscErrorCode DMNetworkGetNumSubNetworks(DM dm,PetscInt *nsubnet,PetscInt *Nsubnet)
61 {
62   DM_Network *network = (DM_Network*)dm->data;
63 
64   PetscFunctionBegin;
65   if (nsubnet) *nsubnet = network->nsubnet;
66   if (Nsubnet) *Nsubnet = network->Nsubnet;
67   PetscFunctionReturn(0);
68 }
69 
70 /*@
71   DMNetworkSetNumSubNetworks - Sets the number of subnetworks
72 
73   Collective on dm
74 
75   Input Parameters:
76 + dm - the dm object
77 . nsubnet - local number of subnetworks
78 - Nsubnet - global number of subnetworks
79 
80    Level: beginner
81 
82 .seealso: DMNetworkCreate(), DMNetworkGetNumSubNetworks()
83 @*/
84 PetscErrorCode DMNetworkSetNumSubNetworks(DM dm,PetscInt nsubnet,PetscInt Nsubnet)
85 {
86   DM_Network     *network = (DM_Network*)dm->data;
87 
88   PetscFunctionBegin;
89   PetscCheck(network->Nsubnet == 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
90 
91   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
92   PetscValidLogicalCollectiveInt(dm,nsubnet,2);
93   PetscValidLogicalCollectiveInt(dm,Nsubnet,3);
94 
95   if (Nsubnet == PETSC_DECIDE) {
96     PetscCheck(nsubnet >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of local subnetworks %" PetscInt_FMT " cannot be less than 0",nsubnet);
97     PetscCall(MPIU_Allreduce(&nsubnet,&Nsubnet,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm)));
98   }
99   PetscCheck(Nsubnet >= 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Number of global subnetworks %" PetscInt_FMT " cannot be less than 1",Nsubnet);
100 
101   network->Nsubnet  = Nsubnet;
102   network->nsubnet  = 0;       /* initia value; will be determind by DMNetworkAddSubnetwork() */
103   PetscCall(PetscCalloc1(Nsubnet,&network->subnet));
104 
105   /* num of shared vertices */
106   network->nsvtx = 0;
107   network->Nsvtx = 0;
108   PetscFunctionReturn(0);
109 }
110 
111 /*@
112   DMNetworkAddSubnetwork - Add a subnetwork
113 
114   Collective on dm
115 
116   Input Parameters:
117 + dm - the dm object
118 . name - name of the subnetwork
119 . ne - number of local edges of this subnetwork
120 - edgelist - list of edges for this subnetwork, this is a one dimensional array with pairs of entries being the two vertices (in global numbering of the vertices) of each edge,
121 $            [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc]
122 
123   Output Parameters:
124 . netnum - global index of the subnetwork
125 
126   Notes:
127   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
128   not be destroyed before the call to DMNetworkLayoutSetUp()
129 
130   A network can comprise of a single subnetwork OR multiple subnetworks. For a single subnetwork, the subnetwork can be read either in serial or parallel. For a multiple subnetworks,
131   each subnetwork topology needs to be set on a unique rank and the communicator size needs to be at least equal to the number of subnetworks.
132 
133   Level: beginner
134 
135   Example usage:
136   Consider the following networks:
137   1) A sigle subnetwork:
138 .vb
139  network 0:
140  rank[0]:
141    v0 --> v2; v1 --> v2
142  rank[1]:
143    v3 --> v5; v4 --> v5
144 .ve
145 
146  The resulting input
147  network 0:
148  rank[0]:
149    ne = 2
150    edgelist = [0 2 | 1 2]
151  rank[1]:
152    ne = 2
153    edgelist = [3 5 | 4 5]
154 
155   2) Two subnetworks:
156 .vb
157  subnetwork 0:
158  rank[0]:
159    v0 --> v2; v2 --> v1; v1 --> v3;
160  subnetwork 1:
161  rank[1]:
162    v0 --> v3; v3 --> v2; v2 --> v1;
163 .ve
164 
165  The resulting input
166  subnetwork 0:
167  rank[0]:
168    ne = 3
169    edgelist = [0 2 | 2 1 | 1 3]
170  rank[1]:
171    ne = 0
172    edgelist = NULL
173 
174  subnetwork 1:
175  rank[0]:
176    ne = 0
177    edgelist = NULL
178  rank[1]:
179    edgelist = [0 3 | 3 2 | 2 1]
180 
181 .seealso: DMNetworkCreate(), DMNetworkSetNumSubnetworks()
182 @*/
183 PetscErrorCode DMNetworkAddSubnetwork(DM dm,const char* name,PetscInt ne,PetscInt edgelist[],PetscInt *netnum)
184 {
185   DM_Network     *network = (DM_Network*)dm->data;
186   PetscInt       i,Nedge,j,Nvtx,nvtx,nvtx_min=-1,nvtx_max=0;
187   PetscBT        table;
188 
189   PetscFunctionBegin;
190   for (i=0; i<ne; i++) {
191     PetscCheck(edgelist[2*i] != edgelist[2*i+1],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Edge %" PetscInt_FMT " has the same vertex %" PetscInt_FMT " at each endpoint",i,edgelist[2*i]);
192   }
193 
194   i = 0;
195   if (ne) nvtx_min = nvtx_max = edgelist[0];
196   for (j=0; j<ne; j++) {
197     nvtx_min = PetscMin(nvtx_min, edgelist[i]);
198     nvtx_max = PetscMax(nvtx_max, edgelist[i]);
199     i++;
200     nvtx_min = PetscMin(nvtx_min, edgelist[i]);
201     nvtx_max = PetscMax(nvtx_max, edgelist[i]);
202     i++;
203   }
204   Nvtx = nvtx_max - nvtx_min + 1; /* approximated total local nvtx for this subnet */
205 
206   /* Get exact local nvtx for this subnet: counting local values between nvtx_min and nvtx_max */
207   PetscCall(PetscBTCreate(Nvtx,&table));
208   PetscCall(PetscBTMemzero(Nvtx,table));
209   i = 0;
210   for (j=0; j<ne; j++) {
211     PetscCall(PetscBTSet(table,edgelist[i++]-nvtx_min));
212     PetscCall(PetscBTSet(table,edgelist[i++]-nvtx_min));
213   }
214   nvtx = 0;
215   for (j=0; j<Nvtx; j++) {
216     if (PetscBTLookup(table,j)) nvtx++;
217   }
218 
219   /* Get global total Nvtx = max(edgelist[])+1 for this subnet */
220   PetscCall(MPIU_Allreduce(&nvtx_max,&Nvtx,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm)));
221   Nvtx++;
222   PetscCall(PetscBTDestroy(&table));
223 
224   /* Get global total Nedge for this subnet */
225   PetscCall(MPIU_Allreduce(&ne,&Nedge,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm)));
226 
227   i = network->nsubnet;
228   if (name) {
229     PetscCall(PetscStrcpy(network->subnet[i].name,name));
230   }
231   network->subnet[i].nvtx     = nvtx; /* include ghost vertices */
232   network->subnet[i].nedge    = ne;
233   network->subnet[i].edgelist = edgelist;
234   network->subnet[i].Nvtx     = Nvtx;
235   network->subnet[i].Nedge    = Nedge;
236 
237   /* ----------------------------------------------------------
238    p=v or e;
239    subnet[0].pStart   = 0
240    subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i])
241    ----------------------------------------------------------------------- */
242   /* GLOBAL subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */
243   network->subnet[i].vStart = network->NVertices;
244   network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx; /* global vEnd of subnet[i] */
245 
246   network->nVertices += nvtx; /* include ghost vertices */
247   network->NVertices += network->subnet[i].Nvtx;
248 
249   /* LOCAL subnet[].eStart and eEnd, used by DMNetworkLayoutSetUp() */
250   network->subnet[i].eStart = network->nEdges;
251   network->subnet[i].eEnd   = network->subnet[i].eStart + ne;
252   network->nEdges += ne;
253   network->NEdges += network->subnet[i].Nedge;
254 
255   PetscCall(PetscStrcpy(network->subnet[i].name,name));
256   if (netnum) *netnum = network->nsubnet;
257   network->nsubnet++;
258   PetscFunctionReturn(0);
259 }
260 
261 /*@C
262   DMNetworkSharedVertexGetInfo - Get info of a shared vertex struct, see petsc/private/dmnetworkimpl.h
263 
264   Not collective
265 
266   Input Parameters:
267 + dm - the DM object
268 - v - vertex point
269 
270   Output Parameters:
271 + gidx - global number of this shared vertex in the internal dmplex
272 . n - number of subnetworks that share this vertex
273 - sv - array of size n: sv[2*i,2*i+1]=(net[i], idx[i]), i=0,...,n-1
274 
275   Level: intermediate
276 
277 .seealso: DMNetworkGetSharedVertices()
278 @*/
279 PetscErrorCode DMNetworkSharedVertexGetInfo(DM dm,PetscInt v,PetscInt *gidx,PetscInt *n,const PetscInt **sv)
280 {
281   DM_Network     *network = (DM_Network*)dm->data;
282   SVtx           *svtx = network->svtx;
283   PetscInt       i,gidx_tmp;
284 
285   PetscFunctionBegin;
286   PetscCall(DMNetworkGetGlobalVertexIndex(dm,v,&gidx_tmp));
287   PetscCall(PetscTableFind(network->svtable,gidx_tmp+1,&i));
288   PetscCheck(i > 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"input vertex is not a shared vertex");
289 
290   i--;
291   if (gidx) *gidx = gidx_tmp;
292   if (n)    *n    = svtx[i].n;
293   if (sv)   *sv   = svtx[i].sv;
294   PetscFunctionReturn(0);
295 }
296 
297 /*
298   VtxGetInfo - Get info of an input vertex=(net,idx)
299 
300   Input Parameters:
301 + Nsvtx - global num of shared vertices
302 . svtx - array of shared vertices (global)
303 - (net,idx) - subnet number and local index for a vertex
304 
305   Output Parameters:
306 + gidx - global index of (net,idx)
307 . svtype - see petsc/private/dmnetworkimpl.h
308 - svtx_idx - ordering in the svtx array
309 */
310 static inline PetscErrorCode VtxGetInfo(PetscInt Nsvtx,SVtx *svtx,PetscInt net,PetscInt idx,PetscInt *gidx,SVtxType *svtype,PetscInt *svtx_idx)
311 {
312   PetscInt i,j,*svto,g_idx;
313   SVtxType vtype;
314 
315   PetscFunctionBegin;
316   if (!Nsvtx) PetscFunctionReturn(0);
317 
318   g_idx = -1;
319   vtype = SVNONE;
320 
321   for (i=0; i<Nsvtx; i++) {
322     if (net == svtx[i].sv[0] && idx == svtx[i].sv[1]) {
323       g_idx = svtx[i].gidx;
324       vtype = SVFROM;
325     } else { /* loop over svtx[i].n */
326       for (j=1; j<svtx[i].n; j++) {
327         svto = svtx[i].sv + 2*j;
328         if (net == svto[0] && idx == svto[1]) {
329           /* input vertex net.idx is a shared to_vertex, output its global index and its svtype */
330           g_idx = svtx[i].gidx; /* output gidx for to_vertex */
331           vtype = SVTO;
332         }
333       }
334     }
335     if (vtype != SVNONE) break;
336   }
337   if (gidx)     *gidx     = g_idx;
338   if (svtype)   *svtype   = vtype;
339   if (svtx_idx) *svtx_idx = i;
340   PetscFunctionReturn(0);
341 }
342 
343 /*
344   TableAddSVtx - Add a new shared vertice from sedgelist[k] to a ctable svta
345 
346   Input:  network, sedgelist, k, svta
347   Output: svta, tdata, ta2sv
348 */
349 static inline PetscErrorCode TableAddSVtx(DM_Network *network,PetscInt *sedgelist,PetscInt k,PetscTable svta,PetscInt* tdata,PetscInt *ta2sv)
350 {
351   PetscInt       net,idx,gidx;
352 
353   PetscFunctionBegin;
354   net = sedgelist[k];
355   idx = sedgelist[k+1];
356   gidx = network->subnet[net].vStart + idx;
357   PetscCall(PetscTableAdd(svta,gidx+1,*tdata+1,INSERT_VALUES));
358 
359   ta2sv[*tdata] = k; /* maps tdata to index of sedgelist */
360   (*tdata)++;
361   PetscFunctionReturn(0);
362 }
363 
364 /*
365   SharedVtxCreate - Create an array of global shared vertices. See SVtx and SVtxType in dmnetworkimpl.h
366 
367   Input:  dm, Nsedgelist, sedgelist
368 
369   Note: Output svtx is organized as
370         sv(net[0],idx[0]) --> sv(net[1],idx[1])
371                           --> sv(net[1],idx[1])
372                           ...
373                           --> sv(net[n-1],idx[n-1])
374         and net[0] < net[1] < ... < net[n-1]
375         where sv[0] has SVFROM type, sv[i], i>0, has SVTO type.
376  */
377 static PetscErrorCode SharedVtxCreate(DM dm,PetscInt Nsedgelist,PetscInt *sedgelist)
378 {
379   SVtx               *svtx = NULL;
380   PetscInt           *sv,k,j,nsv,*tdata,**ta2sv;
381   PetscTable         *svtas;
382   PetscInt           gidx,net,idx,i,nta,ita,idx_from,idx_to,n;
383   DM_Network         *network = (DM_Network*)dm->data;
384   PetscTablePosition ppos;
385 
386   PetscFunctionBegin;
387   /* (1) Crete an array of ctables svtas to map (net,idx) -> gidx; a svtas[] for a shared/merged vertex */
388   PetscCall(PetscCalloc3(Nsedgelist,&svtas,Nsedgelist,&tdata,2*Nsedgelist,&ta2sv));
389 
390   k   = 0;   /* sedgelist vertex counter j = 4*k */
391   nta = 0;   /* num of svta tables created */
392 
393   /* for j=0 */
394   PetscCall(PetscTableCreate(2*Nsedgelist,network->NVertices+1,&svtas[nta]));
395   PetscCall(PetscMalloc1(2*Nsedgelist,&ta2sv[nta]));
396 
397   PetscCall(TableAddSVtx(network,sedgelist,k,svtas[nta],&tdata[nta],ta2sv[nta]));
398   PetscCall(TableAddSVtx(network,sedgelist,k+2,svtas[nta],&tdata[nta],ta2sv[nta]));
399   nta++; k += 4;
400 
401   for (j = 1; j < Nsedgelist; j++) { /* j: sedgelist counter */
402     for (ita = 0; ita < nta; ita++) {
403       /* vfrom */
404       net = sedgelist[k]; idx = sedgelist[k+1];
405       gidx = network->subnet[net].vStart + idx; /* global index of the vertex net.idx before merging shared vertices */
406       PetscCall(PetscTableFind(svtas[ita],gidx+1,&idx_from));
407 
408       /* vto */
409       net = sedgelist[k+2]; idx = sedgelist[k+3];
410       gidx = network->subnet[net].vStart + idx;
411       PetscCall(PetscTableFind(svtas[ita],gidx+1,&idx_to));
412 
413       if (idx_from || idx_to) { /* vfrom or vto is on table svtas[ita] */
414         idx_from--; idx_to--;
415         if (idx_from < 0) { /* vto is on svtas[ita] */
416           PetscCall(TableAddSVtx(network,sedgelist,k,svtas[ita],&tdata[ita],ta2sv[ita]));
417           break;
418         } else if (idx_to < 0) {
419           PetscCall(TableAddSVtx(network,sedgelist,k+2,svtas[ita],&tdata[ita],ta2sv[ita]));
420           break;
421         }
422       }
423     }
424 
425     if (ita == nta) {
426       PetscCall(PetscTableCreate(2*Nsedgelist,network->NVertices+1,&svtas[nta]));
427       PetscCall(PetscMalloc1(2*Nsedgelist, &ta2sv[nta]));
428 
429       PetscCall(TableAddSVtx(network,sedgelist,k,svtas[nta],&tdata[nta],ta2sv[nta]));
430       PetscCall(TableAddSVtx(network,sedgelist,k+2,svtas[nta],&tdata[nta],ta2sv[nta]));
431       nta++;
432     }
433     k += 4;
434   }
435 
436   /* (2) Create svtable for query shared vertices using gidx */
437   PetscCall(PetscTableCreate(nta,network->NVertices+1,&network->svtable));
438 
439   /* (3) Construct svtx from svtas
440      svtx: array of SVtx: sv[0]=(net[0],idx[0]) to vertices sv[k], k=1,...,n-1;
441      net[k], k=0, ...,n-1, are in ascending order */
442   PetscCall(PetscMalloc1(nta,&svtx));
443   for (nsv = 0; nsv < nta; nsv++) {
444     /* for a single svtx, put shared vertices in ascending order of gidx */
445     PetscCall(PetscTableGetCount(svtas[nsv],&n));
446     PetscCall(PetscCalloc1(2*n,&sv));
447     svtx[nsv].sv   = sv;
448     svtx[nsv].n    = n;
449     svtx[nsv].gidx = network->NVertices; /* initialization */
450 
451     PetscCall(PetscTableGetHeadPosition(svtas[nsv],&ppos));
452     for (k=0; k<n; k++) { /* gidx is sorted in ascending order */
453       PetscCall(PetscTableGetNext(svtas[nsv],&ppos,&gidx,&i));
454       gidx--; i--;
455 
456       if (svtx[nsv].gidx > gidx) svtx[nsv].gidx = gidx; /*svtx[nsv].gidx = min(gidx) */
457 
458       j = ta2sv[nsv][i]; /* maps i to index of sedgelist */
459       sv[2*k]   = sedgelist[j];   /* subnet number */
460       sv[2*k+1] = sedgelist[j+1]; /* index on the subnet */
461     }
462 
463     /* Setup svtable for query shared vertices */
464     PetscCall(PetscTableAdd(network->svtable,svtx[nsv].gidx+1,nsv+1,INSERT_VALUES));
465   }
466 
467   for (j=0; j<nta; j++) {
468     PetscCall(PetscTableDestroy(&svtas[j]));
469     PetscCall(PetscFree(ta2sv[j]));
470   }
471   PetscCall(PetscFree3(svtas,tdata,ta2sv));
472 
473   network->Nsvtx = nta;
474   network->svtx  = svtx;
475   PetscFunctionReturn(0);
476 }
477 
478 /*
479   GetEdgelist_Coupling - Get an integrated edgelist for dmplex from user-provided subnet[].edgelist when subnets are coupled by shared vertices
480 
481   Input Parameters:
482 . dm - the dmnetwork object
483 
484    Output Parameters:
485 +  edges - the integrated edgelist for dmplex
486 -  nmerged_ptr - num of vertices being merged
487 */
488 static PetscErrorCode GetEdgelist_Coupling(DM dm,PetscInt *edges,PetscInt *nmerged_ptr)
489 {
490   MPI_Comm       comm;
491   PetscMPIInt    size,rank,*recvcounts=NULL,*displs=NULL;
492   DM_Network     *network = (DM_Network*)dm->data;
493   PetscInt       i,j,ctr,np;
494   PetscInt       *vidxlTog,Nsv,Nsubnet=network->Nsubnet;
495   PetscInt       *sedgelist=network->sedgelist;
496   PetscInt       net,idx,gidx,nmerged,*vrange,gidx_from,net_from,sv_idx;
497   SVtxType       svtype = SVNONE;
498   SVtx           *svtx;
499 
500   PetscFunctionBegin;
501   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
502   PetscCallMPI(MPI_Comm_rank(comm,&rank));
503   PetscCallMPI(MPI_Comm_size(comm,&size));
504 
505   /* (1) Create global svtx[] from sedgelist */
506   /* --------------------------------------- */
507   PetscCall(SharedVtxCreate(dm,network->Nsvtx,sedgelist));
508   Nsv  = network->Nsvtx;
509   svtx = network->svtx;
510 
511   /* (2) Merge shared vto vertices to their vfrom vertex with same global vetex index (gidx) */
512   /* --------------------------------------------------------------------------------------- */
513   /* (2.1) compute vrage[rank]: global index of 1st local vertex in proc[rank] */
514   PetscCall(PetscMalloc4(size+1,&vrange,size,&displs,size,&recvcounts,network->nVertices,&vidxlTog));
515   for (i=0; i<size; i++) {displs[i] = i; recvcounts[i] = 1;}
516 
517   vrange[0] = 0;
518   PetscCallMPI(MPI_Allgatherv(&network->nVertices,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm));
519   for (i=2; i<size+1; i++) vrange[i] += vrange[i-1];
520 
521   /* (2.2) Create vidxlTog: maps UN-MERGED local vertex index i to global index gidx (plex, excluding ghost vertices) */
522   i = 0; gidx = 0;
523   nmerged        = 0; /* local num of merged vertices */
524   network->nsvtx = 0; /* local num of SVtx structs, including ghosts */
525   for (net=0; net<Nsubnet; net++) {
526     for (idx=0; idx<network->subnet[net].Nvtx; idx++) { /* Note: global subnet[net].Nvtx */
527       PetscCall(VtxGetInfo(Nsv,svtx,net,idx,&gidx_from,&svtype,&sv_idx));
528       if (svtype == SVTO) {
529         if (network->subnet[net].nvtx) {/* this proc owns sv_to */
530           net_from = svtx[sv_idx].sv[0]; /* subnet number of its shared vertex */
531           if (network->subnet[net_from].nvtx == 0) {
532             /* this proc does not own v_from, thus a ghost local vertex */
533             network->nsvtx++;
534           }
535           vidxlTog[i++] = gidx_from; /* gidx before merging! Bug??? */
536           nmerged++; /* a shared vertex -- merged */
537         }
538       } else {
539         if (svtype == SVFROM && network->subnet[net].nvtx) {
540           /* this proc owns this v_from, a new local shared vertex */
541           network->nsvtx++;
542         }
543         if (network->subnet[net].nvtx) vidxlTog[i++] = gidx;
544         gidx++;
545       }
546     }
547   }
548 #if defined(PETSC_USE_DEBUG)
549   PetscCheck(i == network->nVertices,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"%" PetscInt_FMT " != %" PetscInt_FMT " nVertices",i,network->nVertices);
550 #endif
551 
552   /* (2.3) Shared vertices in the subnetworks are merged, update global NVertices: np = sum(local nmerged) */
553   PetscCallMPI(MPI_Allreduce(&nmerged,&np,1,MPIU_INT,MPI_SUM,comm));
554   network->NVertices -= np;
555 
556   ctr = 0;
557   for (net=0; net<Nsubnet; net++) {
558     for (j = 0; j < network->subnet[net].nedge; j++) {
559       /* vfrom: */
560       i = network->subnet[net].edgelist[2*j] + (network->subnet[net].vStart - vrange[rank]);
561       edges[2*ctr] = vidxlTog[i];
562 
563       /* vto */
564       i = network->subnet[net].edgelist[2*j+1] + (network->subnet[net].vStart - vrange[rank]);
565       edges[2*ctr+1] = vidxlTog[i];
566       ctr++;
567     }
568   }
569   PetscCall(PetscFree4(vrange,displs,recvcounts,vidxlTog));
570   PetscCall(PetscFree(sedgelist)); /* created in DMNetworkAddSharedVertices() */
571 
572   *nmerged_ptr = nmerged;
573   PetscFunctionReturn(0);
574 }
575 
576 /*@
577   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
578 
579   Not Collective
580 
581   Input Parameters:
582 . dm - the dmnetwork object
583 
584   Notes:
585   This routine should be called after the network sizes and edgelists have been provided. It creates
586   the bare layout of the network and sets up the network to begin insertion of components.
587 
588   All the components should be registered before calling this routine.
589 
590   Level: beginner
591 
592 .seealso: DMNetworkSetNumSubNetworks(), DMNetworkAddSubnetwork()
593 @*/
594 PetscErrorCode DMNetworkLayoutSetUp(DM dm)
595 {
596   DM_Network     *network = (DM_Network*)dm->data;
597   PetscInt       i,j,ctr,Nsubnet=network->Nsubnet,*eowners,np,*edges,*subnetvtx,*subnetedge,e,v,vfrom,vto,net;
598   const PetscInt *cone;
599   MPI_Comm       comm;
600   PetscMPIInt    size,rank;
601   PetscSection   sectiong;
602   PetscInt       nmerged=0;
603 
604   PetscFunctionBegin;
605   PetscCheck(network->nsubnet == Nsubnet,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Must call DMNetworkAddSubnetwork() %" PetscInt_FMT " times",Nsubnet);
606 
607   /* This implementation requires user input each subnet by a single processor when Nsubnet>1, thus subnet[net].nvtx=subnet[net].Nvtx when net>0 */
608   for (net=1; net<Nsubnet; net++) {
609     if (network->subnet[net].nvtx) PetscCheck(network->subnet[net].nvtx == network->subnet[net].Nvtx,PETSC_COMM_SELF,PETSC_ERR_SUP,"subnetwork %" PetscInt_FMT " local num of vertices %" PetscInt_FMT " != %" PetscInt_FMT " global num",net,network->subnet[net].nvtx,network->subnet[net].Nvtx);
610   }
611 
612   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
613   PetscCallMPI(MPI_Comm_rank(comm,&rank));
614   PetscCallMPI(MPI_Comm_size(comm,&size));
615 
616   /* Create LOCAL edgelist in global vertex ordering for the network by concatenating local input edgelists of the subnetworks */
617   PetscCall(PetscCalloc2(2*network->nEdges,&edges,size+1,&eowners));
618 
619   if (network->Nsvtx) { /* subnetworks are coupled via shared vertices */
620     PetscCall(GetEdgelist_Coupling(dm,edges,&nmerged));
621   } else { /* subnetworks are not coupled */
622     /* Create a 0-size svtable for query shared vertices */
623     PetscCall(PetscTableCreate(0,network->NVertices+1,&network->svtable));
624     ctr = 0;
625     for (i=0; i < Nsubnet; i++) {
626       for (j = 0; j < network->subnet[i].nedge; j++) {
627         edges[2*ctr]   = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
628         edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
629         ctr++;
630       }
631     }
632   }
633 
634   /* Create network->plex; One dimensional network, numCorners=2 */
635   PetscCall(DMCreate(comm,&network->plex));
636   PetscCall(DMSetType(network->plex,DMPLEX));
637   PetscCall(DMSetDimension(network->plex,1));
638 
639   if (size == 1) {
640     PetscCall(DMPlexBuildFromCellList(network->plex,network->nEdges,PETSC_DECIDE,2,edges));
641   } else {
642     PetscCall(DMPlexBuildFromCellListParallel(network->plex,network->nEdges,PETSC_DECIDE,PETSC_DECIDE,2,edges,NULL, NULL));
643   }
644 
645   PetscCall(DMPlexGetChart(network->plex,&network->pStart,&network->pEnd));
646   PetscCall(DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd));
647   PetscCall(DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd));
648 
649   PetscCall(PetscSectionCreate(comm,&network->DataSection));
650   PetscCall(PetscSectionCreate(comm,&network->DofSection));
651   PetscCall(PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd));
652   PetscCall(PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd));
653 
654   np = network->pEnd - network->pStart;
655   PetscCall(PetscCalloc2(np,&network->header,np,&network->cvalue));
656   for (i=0; i < np; i++) {
657     network->header[i].maxcomps = 1;
658     PetscCall(SetUpNetworkHeaderComponentValue(dm,&network->header[i],&network->cvalue[i]));
659   }
660 
661   /* Create edge and vertex arrays for the subnetworks
662      This implementation assumes that DMNetwork reads
663      (1) a single subnetwork in parallel; or
664      (2) n subnetworks using n processors, one subnetwork/processor.
665   */
666   PetscCall(PetscCalloc2(network->nEdges,&subnetedge,network->nVertices+network->nsvtx,&subnetvtx)); /* Maps local edge/vertex to local subnetwork's edge/vertex */
667   network->subnetedge = subnetedge;
668   network->subnetvtx  = subnetvtx;
669   for (j=0; j < Nsubnet; j++) {
670     network->subnet[j].edges = subnetedge;
671     subnetedge              += network->subnet[j].nedge;
672 
673     network->subnet[j].vertices = subnetvtx;
674     subnetvtx                  += network->subnet[j].nvtx;
675   }
676   network->svertices = subnetvtx;
677 
678   /* Get edge ownership */
679   np = network->eEnd - network->eStart;
680   PetscCallMPI(MPI_Allgather(&np,1,MPIU_INT,eowners+1,1,MPIU_INT,comm));
681   eowners[0] = 0;
682   for (i=2; i<=size; i++) eowners[i] += eowners[i-1];
683 
684   /* Setup local edge and vertex arrays for subnetworks */
685   e = 0;
686   for (i=0; i < Nsubnet; i++) {
687     ctr = 0;
688     for (j = 0; j < network->subnet[i].nedge; j++) {
689       /* edge e */
690       network->header[e].index    = e + eowners[rank];   /* Global edge index */
691       network->header[e].subnetid = i;
692       network->subnet[i].edges[j] = e;
693 
694       network->header[e].ndata           = 0;
695       network->header[e].offset[0]       = 0;
696       network->header[e].offsetvarrel[0] = 0;
697       PetscCall(PetscSectionAddDof(network->DataSection,e,network->header[e].hsize));
698 
699       /* connected vertices */
700       PetscCall(DMPlexGetCone(network->plex,e,&cone));
701 
702       /* vertex cone[0] */
703       v = cone[0];
704       network->header[v].index     = edges[2*e];  /* Global vertex index */
705       network->header[v].subnetid  = i;           /* Subnetwork id */
706       if (Nsubnet == 1) {
707         network->subnet[i].vertices[v - network->vStart] = v; /* user's subnet[].idx = petsc's v */
708       } else {
709         vfrom = network->subnet[i].edgelist[2*ctr];     /* =subnet[i].idx, Global index! */
710         network->subnet[i].vertices[vfrom] = v; /* user's subnet[].dix = petsc's v */
711       }
712 
713       /* vertex cone[1] */
714       v = cone[1];
715       network->header[v].index    = edges[2*e+1];   /* Global vertex index */
716       network->header[v].subnetid = i;              /* Subnetwork id */
717       if (Nsubnet == 1) {
718         network->subnet[i].vertices[v - network->vStart] = v; /* user's subnet[].idx = petsc's v */
719       } else {
720         vto = network->subnet[i].edgelist[2*ctr+1];     /* =subnet[i].idx, Global index! */
721         network->subnet[i].vertices[vto] = v; /* user's subnet[].dix = petsc's v */
722       }
723 
724       e++; ctr++;
725     }
726   }
727   PetscCall(PetscFree2(edges,eowners));
728 
729   /* Set local vertex array for the subnetworks */
730   j = 0;
731   for (v = network->vStart; v < network->vEnd; v++) {
732     network->header[v].ndata           = 0;
733     network->header[v].offset[0]       = 0;
734     network->header[v].offsetvarrel[0] = 0;
735     PetscCall(PetscSectionAddDof(network->DataSection,v,network->header[e].hsize));
736 
737     /* local shared vertex */
738     PetscCall(PetscTableFind(network->svtable,network->header[v].index+1,&i));
739     if (i) network->svertices[j++] = v;
740   }
741 
742   /* Create a global section to be used by DMNetworkIsGhostVertex() which is a non-collective routine */
743   /* see snes_tutorials_network-ex1_4 */
744   PetscCall(DMGetGlobalSection(network->plex,&sectiong));
745   PetscFunctionReturn(0);
746 }
747 
748 /*@C
749   DMNetworkGetSubnetwork - Returns the information about a requested subnetwork
750 
751   Not collective
752 
753   Input Parameters:
754 + dm - the DM object
755 - netnum - the global index of the subnetwork
756 
757   Output Parameters:
758 + nv - number of vertices (local)
759 . ne - number of edges (local)
760 . vtx - local vertices of the subnetwork
761 - edge - local edges of the subnetwork
762 
763   Notes:
764     Cannot call this routine before DMNetworkLayoutSetup()
765 
766     The local vertices returned on each rank are determined by DMNetwork. The user does not have any control over what vertices are local.
767 
768   Level: intermediate
769 
770 .seealso: DMNetworkCreate(), DMNetworkAddSubnetwork(), DMNetworkLayoutSetUp()
771 @*/
772 PetscErrorCode DMNetworkGetSubnetwork(DM dm,PetscInt netnum,PetscInt *nv,PetscInt *ne,const PetscInt **vtx,const PetscInt **edge)
773 {
774   DM_Network *network = (DM_Network*)dm->data;
775 
776   PetscFunctionBegin;
777   PetscCheck(netnum < network->Nsubnet,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Subnet index %" PetscInt_FMT " exceeds the num of subnets %" PetscInt_FMT "",netnum,network->Nsubnet);
778   if (nv) *nv     = network->subnet[netnum].nvtx;
779   if (ne) *ne     = network->subnet[netnum].nedge;
780   if (vtx) *vtx   = network->subnet[netnum].vertices;
781   if (edge) *edge = network->subnet[netnum].edges;
782   PetscFunctionReturn(0);
783 }
784 
785 /*@
786   DMNetworkAddSharedVertices - Add shared vertices that connect two given subnetworks
787 
788   Collective on dm
789 
790   Input Parameters:
791 + dm - the dm object
792 . anetnum - first subnetwork global numbering returned by DMNetworkAddSubnetwork()
793 . bnetnum - second subnetwork global numbering returned by DMNetworkAddSubnetwork()
794 . nsvtx - number of vertices that are shared by the two subnetworks
795 . asvtx - vertex index in the first subnetwork
796 - bsvtx - vertex index in the second subnetwork
797 
798   Level: beginner
799 
800 .seealso: DMNetworkCreate(), DMNetworkAddSubnetwork(), DMNetworkGetSharedVertices()
801 @*/
802 PetscErrorCode DMNetworkAddSharedVertices(DM dm,PetscInt anetnum,PetscInt bnetnum,PetscInt nsvtx,PetscInt asvtx[],PetscInt bsvtx[])
803 {
804   DM_Network     *network = (DM_Network*)dm->data;
805   PetscInt       i,nsubnet = network->Nsubnet,*sedgelist,Nsvtx=network->Nsvtx;
806 
807   PetscFunctionBegin;
808   PetscCheck(anetnum != bnetnum,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Subnetworks must have different netnum");
809   PetscCheck(anetnum >= 0 && bnetnum >= 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"netnum cannot be negative");
810   if (!Nsvtx) {
811     /* allocate network->sedgelist to hold at most 2*nsubnet pairs of shared vertices */
812     PetscCall(PetscMalloc1(2*4*nsubnet,&network->sedgelist));
813   }
814 
815   sedgelist = network->sedgelist;
816   for (i=0; i<nsvtx; i++) {
817     sedgelist[4*Nsvtx]   = anetnum; sedgelist[4*Nsvtx+1] = asvtx[i];
818     sedgelist[4*Nsvtx+2] = bnetnum; sedgelist[4*Nsvtx+3] = bsvtx[i];
819     Nsvtx++;
820   }
821   PetscCheck(Nsvtx <= 2*nsubnet,PETSC_COMM_SELF,PETSC_ERR_SUP,"allocate more space for coupling edgelist");
822   network->Nsvtx = Nsvtx;
823   PetscFunctionReturn(0);
824 }
825 
826 /*@C
827   DMNetworkGetSharedVertices - Returns the info for the shared vertices
828 
829   Not collective
830 
831   Input Parameter:
832 . dm - the DM object
833 
834   Output Parameters:
835 + nsv - number of local shared vertices
836 - svtx - local shared vertices
837 
838   Notes:
839   Cannot call this routine before DMNetworkLayoutSetup()
840 
841   Level: intermediate
842 
843 .seealso: DMNetworkGetSubnetwork(), DMNetworkLayoutSetUp(), DMNetworkAddSharedVertices()
844 @*/
845 PetscErrorCode DMNetworkGetSharedVertices(DM dm,PetscInt *nsv,const PetscInt **svtx)
846 {
847   DM_Network *net = (DM_Network*)dm->data;
848 
849   PetscFunctionBegin;
850   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
851   if (nsv)  *nsv  = net->nsvtx;
852   if (svtx) *svtx = net->svertices;
853   PetscFunctionReturn(0);
854 }
855 
856 /*@C
857   DMNetworkRegisterComponent - Registers the network component
858 
859   Logically collective on dm
860 
861   Input Parameters:
862 + dm - the network object
863 . name - the component name
864 - size - the storage size in bytes for this component data
865 
866    Output Parameters:
867 .  key - an integer key that defines the component
868 
869    Notes
870    This routine should be called by all processors before calling DMNetworkLayoutSetup().
871 
872    Level: beginner
873 
874 .seealso: DMNetworkCreate(), DMNetworkLayoutSetUp()
875 @*/
876 PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,size_t size,PetscInt *key)
877 {
878   DM_Network            *network = (DM_Network*) dm->data;
879   DMNetworkComponent    *component=NULL,*newcomponent=NULL;
880   PetscBool             flg=PETSC_FALSE;
881   PetscInt              i;
882 
883   PetscFunctionBegin;
884   if (!network->component) {
885     PetscCall(PetscCalloc1(network->max_comps_registered,&network->component));
886   }
887 
888   for (i=0; i < network->ncomponent; i++) {
889     PetscCall(PetscStrcmp(network->component[i].name,name,&flg));
890     if (flg) {
891       *key = i;
892       PetscFunctionReturn(0);
893     }
894   }
895 
896   if (network->ncomponent == network->max_comps_registered) {
897     /* Reached max allowed so resize component */
898     network->max_comps_registered += 2;
899     PetscCall(PetscCalloc1(network->max_comps_registered,&newcomponent));
900     /* Copy over the previous component info */
901     for (i=0; i < network->ncomponent; i++) {
902       PetscCall(PetscStrcpy(newcomponent[i].name,network->component[i].name));
903       newcomponent[i].size = network->component[i].size;
904     }
905     /* Free old one */
906     PetscCall(PetscFree(network->component));
907     /* Update pointer */
908     network->component = newcomponent;
909   }
910 
911   component = &network->component[network->ncomponent];
912 
913   PetscCall(PetscStrcpy(component->name,name));
914   component->size = size/sizeof(DMNetworkComponentGenericDataType);
915   *key = network->ncomponent;
916   network->ncomponent++;
917   PetscFunctionReturn(0);
918 }
919 
920 /*@
921   DMNetworkGetVertexRange - Get the bounds [start, end) for the local vertices
922 
923   Not Collective
924 
925   Input Parameter:
926 . dm - the DMNetwork object
927 
928   Output Parameters:
929 + vStart - the first vertex point
930 - vEnd - one beyond the last vertex point
931 
932   Level: beginner
933 
934 .seealso: DMNetworkGetEdgeRange()
935 @*/
936 PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
937 {
938   DM_Network *network = (DM_Network*)dm->data;
939 
940   PetscFunctionBegin;
941   if (vStart) *vStart = network->vStart;
942   if (vEnd) *vEnd = network->vEnd;
943   PetscFunctionReturn(0);
944 }
945 
946 /*@
947   DMNetworkGetEdgeRange - Get the bounds [start, end) for the local edges
948 
949   Not Collective
950 
951   Input Parameter:
952 . dm - the DMNetwork object
953 
954   Output Parameters:
955 + eStart - The first edge point
956 - eEnd - One beyond the last edge point
957 
958   Level: beginner
959 
960 .seealso: DMNetworkGetVertexRange()
961 @*/
962 PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
963 {
964   DM_Network *network = (DM_Network*)dm->data;
965 
966   PetscFunctionBegin;
967   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
968   if (eStart) *eStart = network->eStart;
969   if (eEnd)   *eEnd   = network->eEnd;
970   PetscFunctionReturn(0);
971 }
972 
973 static PetscErrorCode DMNetworkGetIndex(DM dm,PetscInt p,PetscInt *index)
974 {
975   DM_Network *network = (DM_Network*)dm->data;
976 
977   PetscFunctionBegin;
978   if (network->header) {
979     *index = network->header[p].index;
980   } else {
981     PetscInt                 offsetp;
982     DMNetworkComponentHeader header;
983 
984     PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetp));
985     header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
986     *index = header->index;
987   }
988   PetscFunctionReturn(0);
989 }
990 
991 /*@
992   DMNetworkGetGlobalEdgeIndex - Get the global numbering for the edge on the network
993 
994   Not Collective
995 
996   Input Parameters:
997 + dm - DMNetwork object
998 - p - edge point
999 
1000   Output Parameters:
1001 . index - the global numbering for the edge
1002 
1003   Level: intermediate
1004 
1005 .seealso: DMNetworkGetGlobalVertexIndex()
1006 @*/
1007 PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
1008 {
1009   PetscFunctionBegin;
1010   PetscCall(DMNetworkGetIndex(dm,p,index));
1011   PetscFunctionReturn(0);
1012 }
1013 
1014 /*@
1015   DMNetworkGetGlobalVertexIndex - Get the global numbering for the vertex on the network
1016 
1017   Not Collective
1018 
1019   Input Parameters:
1020 + dm - DMNetwork object
1021 - p  - vertex point
1022 
1023   Output Parameters:
1024 . index - the global numbering for the vertex
1025 
1026   Level: intermediate
1027 
1028 .seealso: DMNetworkGetGlobalEdgeIndex(), DMNetworkGetLocalVertexIndex()
1029 @*/
1030 PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
1031 {
1032   PetscFunctionBegin;
1033   PetscCall(DMNetworkGetIndex(dm,p,index));
1034   PetscFunctionReturn(0);
1035 }
1036 
1037 /*@
1038   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
1039 
1040   Not Collective
1041 
1042   Input Parameters:
1043 + dm - the DMNetwork object
1044 - p - vertex/edge point
1045 
1046   Output Parameters:
1047 . numcomponents - Number of components at the vertex/edge
1048 
1049   Level: beginner
1050 
1051 .seealso: DMNetworkRegisterComponent(), DMNetworkAddComponent()
1052 @*/
1053 PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
1054 {
1055   PetscInt       offset;
1056   DM_Network     *network = (DM_Network*)dm->data;
1057 
1058   PetscFunctionBegin;
1059   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
1060   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
1061   PetscFunctionReturn(0);
1062 }
1063 
1064 /*@
1065   DMNetworkGetLocalVecOffset - Get the offset for accessing the variables associated with a component at the given vertex/edge from the local vector
1066 
1067   Not Collective
1068 
1069   Input Parameters:
1070 + dm - the DMNetwork object
1071 . p - the edge or vertex point
1072 - compnum - component number; use ALL_COMPONENTS if no specific component is requested
1073 
1074   Output Parameters:
1075 . offset - the local offset
1076 
1077   Notes:
1078     These offsets can be passed to MatSetValuesLocal() for matrices obtained with DMCreateMatrix().
1079 
1080     For vectors obtained with DMCreateLocalVector() the offsets can be used with VecSetValues().
1081 
1082     For vectors obtained with DMCreateLocalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1083     the vector values with array[offset].
1084 
1085     For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValuesLocal().
1086 
1087   Level: intermediate
1088 
1089 .seealso: DMGetLocalVector(), DMNetworkGetComponent(), DMNetworkGetGlobalVecOffset(), DMCreateGlobalVector(), VecGetArray(), VecSetValuesLocal(), MatSetValuesLocal()
1090 @*/
1091 PetscErrorCode DMNetworkGetLocalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offset)
1092 {
1093   DM_Network     *network = (DM_Network*)dm->data;
1094   PetscInt       offsetp,offsetd;
1095   DMNetworkComponentHeader header;
1096 
1097   PetscFunctionBegin;
1098   PetscCall(PetscSectionGetOffset(network->plex->localSection,p,&offsetp));
1099   if (compnum == ALL_COMPONENTS) {
1100     *offset = offsetp;
1101     PetscFunctionReturn(0);
1102   }
1103 
1104   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetd));
1105   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd);
1106   *offset = offsetp + header->offsetvarrel[compnum];
1107   PetscFunctionReturn(0);
1108 }
1109 
1110 /*@
1111   DMNetworkGetGlobalVecOffset - Get the global offset for accessing the variables associated with a component for the given vertex/edge from the global vector
1112 
1113   Not Collective
1114 
1115   Input Parameters:
1116 + dm - the DMNetwork object
1117 . p - the edge or vertex point
1118 - compnum - component number; use ALL_COMPONENTS if no specific component is requested
1119 
1120   Output Parameters:
1121 . offsetg - the global offset
1122 
1123   Notes:
1124     These offsets can be passed to MatSetValues() for matrices obtained with DMCreateMatrix().
1125 
1126     For vectors obtained with DMCreateGlobalVector() the offsets can be used with VecSetValues().
1127 
1128     For vectors obtained with DMCreateGlobalVector() and the array obtained with VecGetArray(vec,&array) you can access or set
1129     the vector values with array[offset - rstart] where restart is obtained with VecGetOwnershipRange(v,&rstart,NULL);
1130 
1131   Level: intermediate
1132 
1133 .seealso: DMNetworkGetLocalVecOffset(), DMGetGlobalVector(), DMNetworkGetComponent(), DMCreateGlobalVector(), VecGetArray(), VecSetValues(), MatSetValues()
1134 @*/
1135 PetscErrorCode DMNetworkGetGlobalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offsetg)
1136 {
1137   DM_Network               *network = (DM_Network*)dm->data;
1138   PetscInt                 offsetp,offsetd;
1139   DMNetworkComponentHeader header;
1140 
1141   PetscFunctionBegin;
1142   PetscCall(PetscSectionGetOffset(network->plex->globalSection,p,&offsetp));
1143   if (offsetp < 0) offsetp = -(offsetp + 1); /* Convert to actual global offset for ghost vertex */
1144 
1145   if (compnum == ALL_COMPONENTS) {
1146     *offsetg = offsetp;
1147     PetscFunctionReturn(0);
1148   }
1149   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetd));
1150   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd);
1151   *offsetg = offsetp + header->offsetvarrel[compnum];
1152   PetscFunctionReturn(0);
1153 }
1154 
1155 /*@
1156   DMNetworkGetEdgeOffset - Get the offset for accessing the variables associated with the given edge from the local subvector
1157 
1158   Not Collective
1159 
1160   Input Parameters:
1161 + dm - the DMNetwork object
1162 - p - the edge point
1163 
1164   Output Parameters:
1165 . offset - the offset
1166 
1167   Level: intermediate
1168 
1169 .seealso: DMNetworkGetLocalVecOffset(), DMGetLocalVector()
1170 @*/
1171 PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
1172 {
1173   DM_Network     *network = (DM_Network*)dm->data;
1174 
1175   PetscFunctionBegin;
1176   PetscCall(PetscSectionGetOffset(network->edge.DofSection,p,offset));
1177   PetscFunctionReturn(0);
1178 }
1179 
1180 /*@
1181   DMNetworkGetVertexOffset - Get the offset for accessing the variables associated with the given vertex from the local subvector
1182 
1183   Not Collective
1184 
1185   Input Parameters:
1186 + dm - the DMNetwork object
1187 - p - the vertex point
1188 
1189   Output Parameters:
1190 . offset - the offset
1191 
1192   Level: intermediate
1193 
1194 .seealso: DMNetworkGetEdgeOffset(), DMGetLocalVector()
1195 @*/
1196 PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
1197 {
1198   DM_Network     *network = (DM_Network*)dm->data;
1199 
1200   PetscFunctionBegin;
1201   p -= network->vStart;
1202   PetscCall(PetscSectionGetOffset(network->vertex.DofSection,p,offset));
1203   PetscFunctionReturn(0);
1204 }
1205 
1206 /*@
1207   DMNetworkAddComponent - Adds a network component and number of variables at the given point (vertex/edge)
1208 
1209   Collective on dm
1210 
1211   Input Parameters:
1212 + dm - the DMNetwork
1213 . p - the vertex/edge point. These points are local indices provided by DMNetworkGetSubnetwork()
1214 . componentkey - component key returned while registering the component with DMNetworkRegisterComponent()
1215 . compvalue - pointer to the data structure for the component, or NULL if the component does not require data, this data is not copied so you cannot
1216               free this space until after DMSetUp() is called.
1217 - nvar - number of variables for the component at the vertex/edge point, zero if the component does not introduce any degrees of freedom at the point
1218 
1219   Notes:
1220     The owning rank and any other ranks that have this point as a ghost location must call this routine to add a component and number of variables in the same order at the given point.
1221 
1222     DMNetworkLayoutSetUp() must be called before this routine.
1223 
1224   Developer Notes:
1225      The requirement that all the ranks with access to a vertex (as owner or as ghost) add all the components comes from a limitation of the underlying implementation based on DMPLEX.
1226   Level: beginner
1227 
1228 .seealso: DMNetworkGetComponent(), DMNetworkGetSubnetwork(), DMNetworkIsGhostVertex(), DMNetworkLayoutSetUp()
1229 @*/
1230 PetscErrorCode DMNetworkAddComponent(DM dm,PetscInt p,PetscInt componentkey,void* compvalue,PetscInt nvar)
1231 {
1232   DM_Network               *network = (DM_Network*)dm->data;
1233   DMNetworkComponent       *component = &network->component[componentkey];
1234   DMNetworkComponentHeader header;
1235   DMNetworkComponentValue  cvalue;
1236   PetscInt                 compnum;
1237   PetscInt                 *compsize,*compkey,*compoffset,*compnvar,*compoffsetvarrel;
1238   void*                    *compdata;
1239 
1240   PetscFunctionBegin;
1241   PetscCheck(componentkey >= 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"componentkey %" PetscInt_FMT " cannot be negative. Input a component key returned while registering the component with DMNetworkRegisterComponent()",componentkey);
1242 
1243   /* The owning rank and all ghost ranks add nvar */
1244   PetscCall(PetscSectionAddDof(network->DofSection,p,nvar));
1245 
1246   /* The owning rank and all ghost ranks add a component, including compvalue=NULL */
1247   header = &network->header[p];
1248   cvalue = &network->cvalue[p];
1249   if (header->ndata == header->maxcomps) {
1250     PetscInt additional_size;
1251 
1252     /* Reached limit so resize header component arrays */
1253     header->maxcomps += 2;
1254 
1255     /* Allocate arrays for component information and value */
1256     PetscCall(PetscCalloc5(header->maxcomps,&compsize,header->maxcomps,&compkey,header->maxcomps,&compoffset,header->maxcomps,&compnvar,header->maxcomps,&compoffsetvarrel));
1257     PetscCall(PetscMalloc1(header->maxcomps,&compdata));
1258 
1259     /* Recalculate header size */
1260     header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt);
1261 
1262     header->hsize /= sizeof(DMNetworkComponentGenericDataType);
1263 
1264     /* Copy over component info */
1265     PetscCall(PetscMemcpy(compsize,header->size,header->ndata*sizeof(PetscInt)));
1266     PetscCall(PetscMemcpy(compkey,header->key,header->ndata*sizeof(PetscInt)));
1267     PetscCall(PetscMemcpy(compoffset,header->offset,header->ndata*sizeof(PetscInt)));
1268     PetscCall(PetscMemcpy(compnvar,header->nvar,header->ndata*sizeof(PetscInt)));
1269     PetscCall(PetscMemcpy(compoffsetvarrel,header->offsetvarrel,header->ndata*sizeof(PetscInt)));
1270 
1271     /* Copy over component data pointers */
1272     PetscCall(PetscMemcpy(compdata,cvalue->data,header->ndata*sizeof(void*)));
1273 
1274     /* Free old arrays */
1275     PetscCall(PetscFree5(header->size,header->key,header->offset,header->nvar,header->offsetvarrel));
1276     PetscCall(PetscFree(cvalue->data));
1277 
1278     /* Update pointers */
1279     header->size         = compsize;
1280     header->key          = compkey;
1281     header->offset       = compoffset;
1282     header->nvar         = compnvar;
1283     header->offsetvarrel = compoffsetvarrel;
1284 
1285     cvalue->data = compdata;
1286 
1287     /* Update DataSection Dofs */
1288     /* The dofs for datasection point p equals sizeof the header (i.e. header->hsize) + sizes of the components added at point p. With the resizing of the header, we need to update the dofs for point p. Hence, we add the extra size added for the header */
1289     additional_size = (5*(header->maxcomps - header->ndata)*sizeof(PetscInt))/sizeof(DMNetworkComponentGenericDataType);
1290     PetscCall(PetscSectionAddDof(network->DataSection,p,additional_size));
1291   }
1292   header = &network->header[p];
1293   cvalue = &network->cvalue[p];
1294 
1295   compnum = header->ndata;
1296 
1297   header->size[compnum] = component->size;
1298   PetscCall(PetscSectionAddDof(network->DataSection,p,component->size));
1299   header->key[compnum] = componentkey;
1300   if (compnum != 0) header->offset[compnum] = header->offset[compnum-1] + header->size[compnum-1];
1301   cvalue->data[compnum] = (void*)compvalue;
1302 
1303   /* variables */
1304   header->nvar[compnum] += nvar;
1305   if (compnum != 0) header->offsetvarrel[compnum] = header->offsetvarrel[compnum-1] + header->nvar[compnum-1];
1306 
1307   header->ndata++;
1308   PetscFunctionReturn(0);
1309 }
1310 
1311 /*@
1312   DMNetworkGetComponent - Gets the component key, the component data, and the number of variables at a given network point
1313 
1314   Not Collective
1315 
1316   Input Parameters:
1317 + dm - the DMNetwork object
1318 . p - vertex/edge point
1319 - compnum - component number; use ALL_COMPONENTS if sum up all the components
1320 
1321   Output Parameters:
1322 + compkey - the key obtained when registering the component (use NULL if not required)
1323 . component - the component data (use NULL if not required)
1324 - nvar  - number of variables (use NULL if not required)
1325 
1326   Level: beginner
1327 
1328 .seealso: DMNetworkAddComponent(), DMNetworkGetNumComponents()
1329 @*/
1330 PetscErrorCode DMNetworkGetComponent(DM dm,PetscInt p,PetscInt compnum,PetscInt *compkey,void **component,PetscInt *nvar)
1331 {
1332   DM_Network     *network = (DM_Network*)dm->data;
1333   PetscInt       offset = 0;
1334   DMNetworkComponentHeader header;
1335 
1336   PetscFunctionBegin;
1337   if (compnum == ALL_COMPONENTS) {
1338     PetscCall(PetscSectionGetDof(network->DofSection,p,nvar));
1339     PetscFunctionReturn(0);
1340   }
1341 
1342   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
1343   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
1344 
1345   if (compnum >= 0) {
1346     if (compkey) *compkey = header->key[compnum];
1347     if (component) {
1348       offset += header->hsize+header->offset[compnum];
1349       *component = network->componentdataarray+offset;
1350     }
1351   }
1352 
1353   if (nvar) *nvar = header->nvar[compnum];
1354 
1355   PetscFunctionReturn(0);
1356 }
1357 
1358 /*
1359  Sets up the array that holds the data for all components and its associated section.
1360  It copies the data for all components in a contiguous array called componentdataarray. The component data is stored pointwise with an additional header (metadata) stored for each point. The header has metadata information such as number of components at each point, number of variables for each component, offsets for the components data, etc.
1361 */
1362 PetscErrorCode DMNetworkComponentSetUp(DM dm)
1363 {
1364   DM_Network               *network = (DM_Network*)dm->data;
1365   PetscInt                 arr_size,p,offset,offsetp,ncomp,i,*headerarr;
1366   DMNetworkComponentHeader header;
1367   DMNetworkComponentValue  cvalue;
1368   DMNetworkComponentHeader headerinfo;
1369   DMNetworkComponentGenericDataType *componentdataarray;
1370 
1371   PetscFunctionBegin;
1372   PetscCall(PetscSectionSetUp(network->DataSection));
1373   PetscCall(PetscSectionGetStorageSize(network->DataSection,&arr_size));
1374   /* arr_size+1 fixes pipeline test of opensolaris-misc for src/dm/tests/ex10.c -- Do not know why */
1375   PetscCall(PetscCalloc1(arr_size+1,&network->componentdataarray));
1376 
1377   componentdataarray = network->componentdataarray;
1378   for (p = network->pStart; p < network->pEnd; p++) {
1379     PetscCall(PetscSectionGetOffset(network->DataSection,p,&offsetp));
1380     /* Copy header */
1381     header = &network->header[p];
1382     headerinfo = (DMNetworkComponentHeader)(componentdataarray+offsetp);
1383     PetscCall(PetscMemcpy(headerinfo,header,sizeof(struct _p_DMNetworkComponentHeader)));
1384     headerarr = (PetscInt*)(headerinfo+1);
1385     PetscCall(PetscMemcpy(headerarr,header->size,header->maxcomps*sizeof(PetscInt)));
1386     headerarr += header->maxcomps;
1387     PetscCall(PetscMemcpy(headerarr,header->key,header->maxcomps*sizeof(PetscInt)));
1388     headerarr += header->maxcomps;
1389     PetscCall(PetscMemcpy(headerarr,header->offset,header->maxcomps*sizeof(PetscInt)));
1390     headerarr += header->maxcomps;
1391     PetscCall(PetscMemcpy(headerarr,header->nvar,header->maxcomps*sizeof(PetscInt)));
1392     headerarr += header->maxcomps;
1393     PetscCall(PetscMemcpy(headerarr,header->offsetvarrel,header->maxcomps*sizeof(PetscInt)));
1394 
1395     /* Copy data */
1396     cvalue = &network->cvalue[p];
1397     ncomp  = header->ndata;
1398 
1399     for (i = 0; i < ncomp; i++) {
1400       offset = offsetp + header->hsize + header->offset[i];
1401       PetscCall(PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType)));
1402     }
1403   }
1404   PetscFunctionReturn(0);
1405 }
1406 
1407 /* Sets up the section for dofs. This routine is called during DMSetUp() */
1408 static PetscErrorCode DMNetworkVariablesSetUp(DM dm)
1409 {
1410   DM_Network     *network = (DM_Network*)dm->data;
1411 
1412   PetscFunctionBegin;
1413   PetscCall(PetscSectionSetUp(network->DofSection));
1414   PetscFunctionReturn(0);
1415 }
1416 
1417 /* Get a subsection from a range of points */
1418 static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main,PetscInt pstart,PetscInt pend,PetscSection *subsection)
1419 {
1420   PetscInt       i, nvar;
1421 
1422   PetscFunctionBegin;
1423   PetscCall(PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection));
1424   PetscCall(PetscSectionSetChart(*subsection, 0, pend - pstart));
1425   for (i = pstart; i < pend; i++) {
1426     PetscCall(PetscSectionGetDof(main,i,&nvar));
1427     PetscCall(PetscSectionSetDof(*subsection, i - pstart, nvar));
1428   }
1429 
1430   PetscCall(PetscSectionSetUp(*subsection));
1431   PetscFunctionReturn(0);
1432 }
1433 
1434 /* Create a submap of points with a GlobalToLocal structure */
1435 static PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
1436 {
1437   PetscInt       i, *subpoints;
1438 
1439   PetscFunctionBegin;
1440   /* Create index sets to map from "points" to "subpoints" */
1441   PetscCall(PetscMalloc1(pend - pstart, &subpoints));
1442   for (i = pstart; i < pend; i++) {
1443     subpoints[i - pstart] = i;
1444   }
1445   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map));
1446   PetscCall(PetscFree(subpoints));
1447   PetscFunctionReturn(0);
1448 }
1449 
1450 /*@
1451   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute
1452 
1453   Collective on dm
1454 
1455   Input Parameters:
1456 . dm - the DMNetworkObject
1457 
1458   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
1459 
1460   points = [0 1 2 3 4 5 6]
1461 
1462   where edges = [0,1,2,3] and vertices = [4,5,6]. The new orderings will be specific to the subset (i.e vertices = [0,1,2] <- [4,5,6]).
1463 
1464   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
1465 
1466   Level: intermediate
1467 
1468 @*/
1469 PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
1470 {
1471   MPI_Comm       comm;
1472   PetscMPIInt    size;
1473   DM_Network     *network = (DM_Network*)dm->data;
1474 
1475   PetscFunctionBegin;
1476   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
1477   PetscCallMPI(MPI_Comm_size(comm, &size));
1478 
1479   /* Create maps for vertices and edges */
1480   PetscCall(DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping));
1481   PetscCall(DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping));
1482 
1483   /* Create local sub-sections */
1484   PetscCall(DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection));
1485   PetscCall(DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection));
1486 
1487   if (size > 1) {
1488     PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf));
1489 
1490     PetscCall(PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection));
1491     PetscCall(PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf));
1492     PetscCall(PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection));
1493   } else {
1494     /* create structures for vertex */
1495     PetscCall(PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection));
1496     /* create structures for edge */
1497     PetscCall(PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection));
1498   }
1499 
1500   /* Add viewers */
1501   PetscCall(PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section"));
1502   PetscCall(PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section"));
1503   PetscCall(PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view"));
1504   PetscCall(PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view"));
1505   PetscFunctionReturn(0);
1506 }
1507 
1508 /*
1509    Setup a lookup btable for the input v's owning subnetworks
1510    - add all owing subnetworks that connect to this v to the btable
1511      vertex_subnetid = supportingedge_subnetid
1512 */
1513 static inline PetscErrorCode SetSubnetIdLookupBT(DM dm,PetscInt v,PetscInt Nsubnet,PetscBT btable)
1514 {
1515   PetscInt       e,nedges,offset;
1516   const PetscInt *edges;
1517   DM_Network     *newDMnetwork = (DM_Network*)dm->data;
1518   DMNetworkComponentHeader header;
1519 
1520   PetscFunctionBegin;
1521   PetscCall(PetscBTMemzero(Nsubnet,btable));
1522   PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges));
1523   for (e=0; e<nedges; e++) {
1524     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,edges[e],&offset));
1525     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1526     PetscCall(PetscBTSet(btable,header->subnetid));
1527   }
1528   PetscFunctionReturn(0);
1529 }
1530 
1531 /*@
1532   DMNetworkDistribute - Distributes the network and moves associated component data
1533 
1534   Collective
1535 
1536   Input Parameters:
1537 + DM - the DMNetwork object
1538 - overlap - the overlap of partitions, 0 is the default
1539 
1540   Options Database Key:
1541 + -dmnetwork_view - Calls DMView() at the conclusion of DMSetUp()
1542 - -dmnetwork_view_distributed - Calls DMView() at the conclusion of DMNetworkDistribute()
1543 
1544   Notes:
1545   Distributes the network with <overlap>-overlapping partitioning of the edges.
1546 
1547   Level: intermediate
1548 
1549 .seealso: DMNetworkCreate()
1550 @*/
1551 PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
1552 {
1553   MPI_Comm       comm;
1554   PetscMPIInt    size;
1555   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
1556   DM_Network     *newDMnetwork;
1557   PetscSF        pointsf=NULL;
1558   DM             newDM;
1559   PetscInt       j,e,v,offset,*subnetvtx,*subnetedge,Nsubnet,gidx,svtx_idx,nv;
1560   PetscInt       net,*sv;
1561   PetscBT        btable;
1562   PetscPartitioner         part;
1563   DMNetworkComponentHeader header;
1564 
1565   PetscFunctionBegin;
1566   PetscValidPointer(dm,1);
1567   PetscValidHeaderSpecific(*dm,DM_CLASSID,1);
1568   PetscCall(PetscObjectGetComm((PetscObject)*dm,&comm));
1569   PetscCallMPI(MPI_Comm_size(comm, &size));
1570   if (size == 1) PetscFunctionReturn(0);
1571 
1572   PetscCheck(!overlap,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"overlap %" PetscInt_FMT " != 0 is not supported yet",overlap);
1573 
1574   /* This routine moves the component data to the appropriate processors. It makes use of the DataSection and the componentdataarray to move the component data to appropriate processors and returns a new DataSection and new componentdataarray. */
1575   PetscCall(DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM));
1576   newDMnetwork = (DM_Network*)newDM->data;
1577   newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered;
1578   PetscCall(PetscMalloc1(newDMnetwork->max_comps_registered,&newDMnetwork->component));
1579 
1580   /* Enable runtime options for petscpartitioner */
1581   PetscCall(DMPlexGetPartitioner(oldDMnetwork->plex,&part));
1582   PetscCall(PetscPartitionerSetFromOptions(part));
1583 
1584   /* Distribute plex dm */
1585   PetscCall(DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex));
1586 
1587   /* Distribute dof section */
1588   PetscCall(PetscSectionCreate(comm,&newDMnetwork->DofSection));
1589   PetscCall(PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection));
1590 
1591   /* Distribute data and associated section */
1592   PetscCall(PetscSectionCreate(comm,&newDMnetwork->DataSection));
1593   PetscCall(DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray));
1594 
1595   PetscCall(PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd));
1596   PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd));
1597   PetscCall(DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd));
1598   newDMnetwork->nEdges    = newDMnetwork->eEnd - newDMnetwork->eStart;
1599   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
1600   newDMnetwork->NVertices = oldDMnetwork->NVertices;
1601   newDMnetwork->NEdges    = oldDMnetwork->NEdges;
1602   newDMnetwork->svtable   = oldDMnetwork->svtable; /* global table! */
1603   oldDMnetwork->svtable   = NULL;
1604 
1605   /* Set Dof section as the section for dm */
1606   PetscCall(DMSetLocalSection(newDMnetwork->plex,newDMnetwork->DofSection));
1607   PetscCall(DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection));
1608 
1609   /* Setup subnetwork info in the newDM */
1610   newDMnetwork->Nsubnet = oldDMnetwork->Nsubnet;
1611   newDMnetwork->Nsvtx   = oldDMnetwork->Nsvtx;
1612   oldDMnetwork->Nsvtx   = 0;
1613   newDMnetwork->svtx    = oldDMnetwork->svtx; /* global vertices! */
1614   oldDMnetwork->svtx    = NULL;
1615   PetscCall(PetscCalloc1(newDMnetwork->Nsubnet,&newDMnetwork->subnet));
1616 
1617   /* Copy over the global number of vertices and edges in each subnetwork.
1618      Note: these are calculated in DMNetworkLayoutSetUp()
1619   */
1620   Nsubnet = newDMnetwork->Nsubnet;
1621   for (j = 0; j < Nsubnet; j++) {
1622     newDMnetwork->subnet[j].Nvtx  = oldDMnetwork->subnet[j].Nvtx;
1623     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1624   }
1625 
1626   /* Count local nedges for subnetworks */
1627   for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) {
1628     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset));
1629     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1630 
1631     /* Update pointers */
1632     header->size          = (PetscInt*)(header + 1);
1633     header->key           = header->size   + header->maxcomps;
1634     header->offset        = header->key    + header->maxcomps;
1635     header->nvar          = header->offset + header->maxcomps;
1636     header->offsetvarrel  = header->nvar   + header->maxcomps;
1637 
1638     newDMnetwork->subnet[header->subnetid].nedge++;
1639   }
1640 
1641   /* Setup a btable to keep track subnetworks owned by this process at a shared vertex */
1642   if (newDMnetwork->Nsvtx) {
1643     PetscCall(PetscBTCreate(Nsubnet,&btable));
1644   }
1645 
1646   /* Count local nvtx for subnetworks */
1647   for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) {
1648     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset));
1649     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1650 
1651     /* Update pointers */
1652     header->size          = (PetscInt*)(header + 1);
1653     header->key           = header->size   + header->maxcomps;
1654     header->offset        = header->key    + header->maxcomps;
1655     header->nvar          = header->offset + header->maxcomps;
1656     header->offsetvarrel  = header->nvar   + header->maxcomps;
1657 
1658     /* shared vertices: use gidx=header->index to check if v is a shared vertex */
1659     gidx = header->index;
1660     PetscCall(PetscTableFind(newDMnetwork->svtable,gidx+1,&svtx_idx));
1661     svtx_idx--;
1662 
1663     if (svtx_idx < 0) { /* not a shared vertex */
1664       newDMnetwork->subnet[header->subnetid].nvtx++;
1665     } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */
1666       /* Setup a lookup btable for this v's owning subnetworks */
1667       PetscCall(SetSubnetIdLookupBT(newDM,v,Nsubnet,btable));
1668 
1669       for (j=0; j<newDMnetwork->svtx[svtx_idx].n; j++) {
1670         sv  = newDMnetwork->svtx[svtx_idx].sv + 2*j;
1671         net = sv[0];
1672         if (PetscBTLookup(btable,net)) newDMnetwork->subnet[net].nvtx++; /* sv is on net owned by this proces */
1673       }
1674     }
1675   }
1676 
1677   /* Get total local nvtx for subnetworks */
1678   nv = 0;
1679   for (j=0; j<Nsubnet; j++) nv += newDMnetwork->subnet[j].nvtx;
1680   nv += newDMnetwork->Nsvtx;
1681 
1682   /* Now create the vertices and edge arrays for the subnetworks */
1683   PetscCall(PetscCalloc2(newDMnetwork->nEdges,&subnetedge,nv,&subnetvtx)); /* Maps local vertex to local subnetwork's vertex */
1684   newDMnetwork->subnetedge = subnetedge;
1685   newDMnetwork->subnetvtx  = subnetvtx;
1686   for (j=0; j < newDMnetwork->Nsubnet; j++) {
1687     newDMnetwork->subnet[j].edges = subnetedge;
1688     subnetedge                   += newDMnetwork->subnet[j].nedge;
1689 
1690     newDMnetwork->subnet[j].vertices = subnetvtx;
1691     subnetvtx                       += newDMnetwork->subnet[j].nvtx;
1692 
1693     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. These get updated when the vertices and edges are added. */
1694     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1695   }
1696   newDMnetwork->svertices = subnetvtx;
1697 
1698   /* Set the edges and vertices in each subnetwork */
1699   for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) {
1700     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset));
1701     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1702     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e;
1703   }
1704 
1705   nv = 0;
1706   for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) {
1707     PetscCall(PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset));
1708     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);
1709 
1710     /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */
1711     PetscCall(PetscTableFind(newDMnetwork->svtable,header->index+1,&svtx_idx));
1712     svtx_idx--;
1713     if (svtx_idx < 0) {
1714       newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v;
1715     } else { /* a shared vertex */
1716       newDMnetwork->svertices[nv++] = v;
1717 
1718       /* Setup a lookup btable for this v's owning subnetworks */
1719       PetscCall(SetSubnetIdLookupBT(newDM,v,Nsubnet,btable));
1720 
1721       for (j=0; j<newDMnetwork->svtx[svtx_idx].n; j++) {
1722         sv  = newDMnetwork->svtx[svtx_idx].sv + 2*j;
1723         net = sv[0];
1724         if (PetscBTLookup(btable,net))
1725           newDMnetwork->subnet[net].vertices[newDMnetwork->subnet[net].nvtx++] = v;
1726       }
1727     }
1728   }
1729   newDMnetwork->nsvtx = nv;   /* num of local shared vertices */
1730 
1731   newDM->setupcalled = (*dm)->setupcalled;
1732   newDMnetwork->distributecalled = PETSC_TRUE;
1733 
1734   /* Free spaces */
1735   PetscCall(PetscSFDestroy(&pointsf));
1736   PetscCall(DMDestroy(dm));
1737   if (newDMnetwork->Nsvtx) {
1738     PetscCall(PetscBTDestroy(&btable));
1739   }
1740 
1741   /* View distributed dmnetwork */
1742   PetscCall(DMViewFromOptions(newDM,NULL,"-dmnetwork_view_distributed"));
1743 
1744   *dm  = newDM;
1745   PetscFunctionReturn(0);
1746 }
1747 
1748 /*@C
1749   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering
1750 
1751  Collective
1752 
1753   Input Parameters:
1754 + mainSF - the original SF structure
1755 - map - a ISLocalToGlobal mapping that contains the subset of points
1756 
1757   Output Parameter:
1758 . subSF - a subset of the mainSF for the desired subset.
1759 
1760   Level: intermediate
1761 @*/
1762 PetscErrorCode PetscSFGetSubSF(PetscSF mainsf,ISLocalToGlobalMapping map,PetscSF *subSF)
1763 {
1764   PetscInt              nroots, nleaves, *ilocal_sub;
1765   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
1766   PetscInt              *local_points, *remote_points;
1767   PetscSFNode           *iremote_sub;
1768   const PetscInt        *ilocal;
1769   const PetscSFNode     *iremote;
1770 
1771   PetscFunctionBegin;
1772   PetscCall(PetscSFGetGraph(mainsf,&nroots,&nleaves,&ilocal,&iremote));
1773 
1774   /* Look for leaves that pertain to the subset of points. Get the local ordering */
1775   PetscCall(PetscMalloc1(nleaves,&ilocal_map));
1776   PetscCall(ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map));
1777   for (i = 0; i < nleaves; i++) {
1778     if (ilocal_map[i] != -1) nleaves_sub += 1;
1779   }
1780   /* Re-number ilocal with subset numbering. Need information from roots */
1781   PetscCall(PetscMalloc2(nroots,&local_points,nroots,&remote_points));
1782   for (i = 0; i < nroots; i++) local_points[i] = i;
1783   PetscCall(ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points));
1784   PetscCall(PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE));
1785   PetscCall(PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE));
1786   /* Fill up graph using local (that is, local to the subset) numbering. */
1787   PetscCall(PetscMalloc1(nleaves_sub,&ilocal_sub));
1788   PetscCall(PetscMalloc1(nleaves_sub,&iremote_sub));
1789   nleaves_sub = 0;
1790   for (i = 0; i < nleaves; i++) {
1791     if (ilocal_map[i] != -1) {
1792       ilocal_sub[nleaves_sub] = ilocal_map[i];
1793       iremote_sub[nleaves_sub].rank = iremote[i].rank;
1794       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
1795       nleaves_sub += 1;
1796     }
1797   }
1798   PetscCall(PetscFree2(local_points,remote_points));
1799   PetscCall(ISLocalToGlobalMappingGetSize(map,&nroots_sub));
1800 
1801   /* Create new subSF */
1802   PetscCall(PetscSFCreate(PETSC_COMM_WORLD,subSF));
1803   PetscCall(PetscSFSetFromOptions(*subSF));
1804   PetscCall(PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES));
1805   PetscCall(PetscFree(ilocal_map));
1806   PetscCall(PetscFree(iremote_sub));
1807   PetscFunctionReturn(0);
1808 }
1809 
1810 /*@C
1811   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
1812 
1813   Not Collective
1814 
1815   Input Parameters:
1816 + dm - the DMNetwork object
1817 - p  - the vertex point
1818 
1819   Output Parameters:
1820 + nedges - number of edges connected to this vertex point
1821 - edges  - list of edge points
1822 
1823   Level: beginner
1824 
1825   Fortran Notes:
1826   Since it returns an array, this routine is only available in Fortran 90, and you must
1827   include petsc.h90 in your code.
1828 
1829 .seealso: DMNetworkCreate(), DMNetworkGetConnectedVertices()
1830 @*/
1831 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
1832 {
1833   DM_Network     *network = (DM_Network*)dm->data;
1834 
1835   PetscFunctionBegin;
1836   PetscCall(DMPlexGetSupportSize(network->plex,vertex,nedges));
1837   if (edges) PetscCall(DMPlexGetSupport(network->plex,vertex,edges));
1838   PetscFunctionReturn(0);
1839 }
1840 
1841 /*@C
1842   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
1843 
1844   Not Collective
1845 
1846   Input Parameters:
1847 + dm - the DMNetwork object
1848 - p - the edge point
1849 
1850   Output Parameters:
1851 . vertices - vertices connected to this edge
1852 
1853   Level: beginner
1854 
1855   Fortran Notes:
1856   Since it returns an array, this routine is only available in Fortran 90, and you must
1857   include petsc.h90 in your code.
1858 
1859 .seealso: DMNetworkCreate(), DMNetworkGetSupportingEdges()
1860 @*/
1861 PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
1862 {
1863   DM_Network     *network = (DM_Network*)dm->data;
1864 
1865   PetscFunctionBegin;
1866   PetscCall(DMPlexGetCone(network->plex,edge,vertices));
1867   PetscFunctionReturn(0);
1868 }
1869 
1870 /*@
1871   DMNetworkIsSharedVertex - Returns TRUE if the vertex is shared by subnetworks
1872 
1873   Not Collective
1874 
1875   Input Parameters:
1876 + dm - the DMNetwork object
1877 - p - the vertex point
1878 
1879   Output Parameter:
1880 . flag - TRUE if the vertex is shared by subnetworks
1881 
1882   Level: beginner
1883 
1884 .seealso: DMNetworkAddSharedVertices(), DMNetworkIsGhostVertex()
1885 @*/
1886 PetscErrorCode DMNetworkIsSharedVertex(DM dm,PetscInt p,PetscBool *flag)
1887 {
1888   PetscInt       i;
1889 
1890   PetscFunctionBegin;
1891   *flag = PETSC_FALSE;
1892 
1893   if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */
1894     DM_Network     *network = (DM_Network*)dm->data;
1895     PetscInt       gidx;
1896     PetscCall(DMNetworkGetGlobalVertexIndex(dm,p,&gidx));
1897     PetscCall(PetscTableFind(network->svtable,gidx+1,&i));
1898     if (i) *flag = PETSC_TRUE;
1899   } else { /* would be removed? */
1900     PetscInt       nv;
1901     const PetscInt *vtx;
1902     PetscCall(DMNetworkGetSharedVertices(dm,&nv,&vtx));
1903     for (i=0; i<nv; i++) {
1904       if (p == vtx[i]) {
1905         *flag = PETSC_TRUE;
1906         break;
1907       }
1908     }
1909   }
1910   PetscFunctionReturn(0);
1911 }
1912 
1913 /*@
1914   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
1915 
1916   Not Collective
1917 
1918   Input Parameters:
1919 + dm - the DMNetwork object
1920 - p - the vertex point
1921 
1922   Output Parameter:
1923 . isghost - TRUE if the vertex is a ghost point
1924 
1925   Level: beginner
1926 
1927 .seealso: DMNetworkGetConnectedVertices(), DMNetworkGetVertexRange(), DMNetworkIsSharedVertex()
1928 @*/
1929 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
1930 {
1931   DM_Network     *network = (DM_Network*)dm->data;
1932   PetscInt       offsetg;
1933   PetscSection   sectiong;
1934 
1935   PetscFunctionBegin;
1936   *isghost = PETSC_FALSE;
1937   PetscCall(DMGetGlobalSection(network->plex,&sectiong));
1938   PetscCall(PetscSectionGetOffset(sectiong,p,&offsetg));
1939   if (offsetg < 0) *isghost = PETSC_TRUE;
1940   PetscFunctionReturn(0);
1941 }
1942 
1943 PetscErrorCode DMSetUp_Network(DM dm)
1944 {
1945   DM_Network     *network=(DM_Network*)dm->data;
1946 
1947   PetscFunctionBegin;
1948   PetscCall(DMNetworkComponentSetUp(dm));
1949   PetscCall(DMNetworkVariablesSetUp(dm));
1950 
1951   PetscCall(DMSetLocalSection(network->plex,network->DofSection));
1952   PetscCall(DMGetGlobalSection(network->plex,&network->GlobalDofSection));
1953 
1954   dm->setupcalled = PETSC_TRUE;
1955 
1956   /* View dmnetwork */
1957   PetscCall(DMViewFromOptions(dm,NULL,"-dmnetwork_view"));
1958   PetscFunctionReturn(0);
1959 }
1960 
1961 /*@
1962   DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
1963       -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
1964 
1965   Collective
1966 
1967   Input Parameters:
1968 + dm - the DMNetwork object
1969 . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
1970 - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
1971 
1972  Level: intermediate
1973 
1974 @*/
1975 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
1976 {
1977   DM_Network     *network=(DM_Network*)dm->data;
1978   PetscInt       nVertices = network->nVertices;
1979 
1980   PetscFunctionBegin;
1981   network->userEdgeJacobian   = eflg;
1982   network->userVertexJacobian = vflg;
1983 
1984   if (eflg && !network->Je) {
1985     PetscCall(PetscCalloc1(3*network->nEdges,&network->Je));
1986   }
1987 
1988   if (vflg && !network->Jv && nVertices) {
1989     PetscInt       i,*vptr,nedges,vStart=network->vStart;
1990     PetscInt       nedges_total;
1991     const PetscInt *edges;
1992 
1993     /* count nvertex_total */
1994     nedges_total = 0;
1995     PetscCall(PetscMalloc1(nVertices+1,&vptr));
1996 
1997     vptr[0] = 0;
1998     for (i=0; i<nVertices; i++) {
1999       PetscCall(DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges));
2000       nedges_total += nedges;
2001       vptr[i+1] = vptr[i] + 2*nedges + 1;
2002     }
2003 
2004     PetscCall(PetscCalloc1(2*nedges_total+nVertices,&network->Jv));
2005     network->Jvptr = vptr;
2006   }
2007   PetscFunctionReturn(0);
2008 }
2009 
2010 /*@
2011   DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
2012 
2013   Not Collective
2014 
2015   Input Parameters:
2016 + dm - the DMNetwork object
2017 . p - the edge point
2018 - J - array (size = 3) of Jacobian submatrices for this edge point:
2019         J[0]: this edge
2020         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
2021 
2022   Level: advanced
2023 
2024 .seealso: DMNetworkVertexSetMatrix()
2025 @*/
2026 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
2027 {
2028   DM_Network *network=(DM_Network*)dm->data;
2029 
2030   PetscFunctionBegin;
2031   PetscCheck(network->Je,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
2032 
2033   if (J) {
2034     network->Je[3*p]   = J[0];
2035     network->Je[3*p+1] = J[1];
2036     network->Je[3*p+2] = J[2];
2037   }
2038   PetscFunctionReturn(0);
2039 }
2040 
2041 /*@
2042   DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
2043 
2044   Not Collective
2045 
2046   Input Parameters:
2047 + dm - The DMNetwork object
2048 . p - the vertex point
2049 - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
2050         J[0]:       this vertex
2051         J[1+2*i]:   i-th supporting edge
2052         J[1+2*i+1]: i-th connected vertex
2053 
2054   Level: advanced
2055 
2056 .seealso: DMNetworkEdgeSetMatrix()
2057 @*/
2058 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
2059 {
2060   DM_Network     *network=(DM_Network*)dm->data;
2061   PetscInt       i,*vptr,nedges,vStart=network->vStart;
2062   const PetscInt *edges;
2063 
2064   PetscFunctionBegin;
2065   PetscCheck(network->Jv,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
2066 
2067   if (J) {
2068     vptr = network->Jvptr;
2069     network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
2070 
2071     /* Set Jacobian for each supporting edge and connected vertex */
2072     PetscCall(DMNetworkGetSupportingEdges(dm,p,&nedges,&edges));
2073     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
2074   }
2075   PetscFunctionReturn(0);
2076 }
2077 
2078 static inline PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2079 {
2080   PetscInt       j;
2081   PetscScalar    val=(PetscScalar)ncols;
2082 
2083   PetscFunctionBegin;
2084   if (!ghost) {
2085     for (j=0; j<nrows; j++) {
2086       PetscCall(VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES));
2087     }
2088   } else {
2089     for (j=0; j<nrows; j++) {
2090       PetscCall(VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES));
2091     }
2092   }
2093   PetscFunctionReturn(0);
2094 }
2095 
2096 static inline PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2097 {
2098   PetscInt       j,ncols_u;
2099   PetscScalar    val;
2100 
2101   PetscFunctionBegin;
2102   if (!ghost) {
2103     for (j=0; j<nrows; j++) {
2104       PetscCall(MatGetRow(Ju,j,&ncols_u,NULL,NULL));
2105       val = (PetscScalar)ncols_u;
2106       PetscCall(VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES));
2107       PetscCall(MatRestoreRow(Ju,j,&ncols_u,NULL,NULL));
2108     }
2109   } else {
2110     for (j=0; j<nrows; j++) {
2111       PetscCall(MatGetRow(Ju,j,&ncols_u,NULL,NULL));
2112       val = (PetscScalar)ncols_u;
2113       PetscCall(VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES));
2114       PetscCall(MatRestoreRow(Ju,j,&ncols_u,NULL,NULL));
2115     }
2116   }
2117   PetscFunctionReturn(0);
2118 }
2119 
2120 static inline PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
2121 {
2122   PetscFunctionBegin;
2123   if (Ju) {
2124     PetscCall(MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz));
2125   } else {
2126     PetscCall(MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz));
2127   }
2128   PetscFunctionReturn(0);
2129 }
2130 
2131 static inline PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2132 {
2133   PetscInt       j,*cols;
2134   PetscScalar    *zeros;
2135 
2136   PetscFunctionBegin;
2137   PetscCall(PetscCalloc2(ncols,&cols,nrows*ncols,&zeros));
2138   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
2139   PetscCall(MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES));
2140   PetscCall(PetscFree2(cols,zeros));
2141   PetscFunctionReturn(0);
2142 }
2143 
2144 static inline PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2145 {
2146   PetscInt       j,M,N,row,col,ncols_u;
2147   const PetscInt *cols;
2148   PetscScalar    zero=0.0;
2149 
2150   PetscFunctionBegin;
2151   PetscCall(MatGetSize(Ju,&M,&N));
2152   PetscCheck(nrows == M && ncols == N,PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%" PetscInt_FMT " by %" PetscInt_FMT " must equal %" PetscInt_FMT " by %" PetscInt_FMT "",nrows,ncols,M,N);
2153 
2154   for (row=0; row<nrows; row++) {
2155     PetscCall(MatGetRow(Ju,row,&ncols_u,&cols,NULL));
2156     for (j=0; j<ncols_u; j++) {
2157       col = cols[j] + cstart;
2158       PetscCall(MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES));
2159     }
2160     PetscCall(MatRestoreRow(Ju,row,&ncols_u,&cols,NULL));
2161   }
2162   PetscFunctionReturn(0);
2163 }
2164 
2165 static inline PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
2166 {
2167   PetscFunctionBegin;
2168   if (Ju) {
2169     PetscCall(MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J));
2170   } else {
2171     PetscCall(MatSetDenseblock_private(nrows,rows,ncols,cstart,J));
2172   }
2173   PetscFunctionReturn(0);
2174 }
2175 
2176 /* Creates a GlobalToLocal mapping with a Local and Global section. This is akin to the routine DMGetLocalToGlobalMapping but without the need of providing a dm.
2177 */
2178 PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
2179 {
2180   PetscInt       i,size,dof;
2181   PetscInt       *glob2loc;
2182 
2183   PetscFunctionBegin;
2184   PetscCall(PetscSectionGetStorageSize(localsec,&size));
2185   PetscCall(PetscMalloc1(size,&glob2loc));
2186 
2187   for (i = 0; i < size; i++) {
2188     PetscCall(PetscSectionGetOffset(globalsec,i,&dof));
2189     dof = (dof >= 0) ? dof : -(dof + 1);
2190     glob2loc[i] = dof;
2191   }
2192 
2193   PetscCall(ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog));
2194 #if 0
2195   PetscCall(PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD));
2196 #endif
2197   PetscFunctionReturn(0);
2198 }
2199 
2200 #include <petsc/private/matimpl.h>
2201 
2202 PetscErrorCode DMCreateMatrix_Network_Nest(DM dm,Mat *J)
2203 {
2204   DM_Network     *network = (DM_Network*)dm->data;
2205   PetscInt       eDof,vDof;
2206   Mat            j11,j12,j21,j22,bA[2][2];
2207   MPI_Comm       comm;
2208   ISLocalToGlobalMapping eISMap,vISMap;
2209 
2210   PetscFunctionBegin;
2211   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
2212 
2213   PetscCall(PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof));
2214   PetscCall(PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof));
2215 
2216   PetscCall(MatCreate(comm, &j11));
2217   PetscCall(MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE));
2218   PetscCall(MatSetType(j11, MATMPIAIJ));
2219 
2220   PetscCall(MatCreate(comm, &j12));
2221   PetscCall(MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE));
2222   PetscCall(MatSetType(j12, MATMPIAIJ));
2223 
2224   PetscCall(MatCreate(comm, &j21));
2225   PetscCall(MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE));
2226   PetscCall(MatSetType(j21, MATMPIAIJ));
2227 
2228   PetscCall(MatCreate(comm, &j22));
2229   PetscCall(MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE));
2230   PetscCall(MatSetType(j22, MATMPIAIJ));
2231 
2232   bA[0][0] = j11;
2233   bA[0][1] = j12;
2234   bA[1][0] = j21;
2235   bA[1][1] = j22;
2236 
2237   PetscCall(CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap));
2238   PetscCall(CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap));
2239 
2240   PetscCall(MatSetLocalToGlobalMapping(j11,eISMap,eISMap));
2241   PetscCall(MatSetLocalToGlobalMapping(j12,eISMap,vISMap));
2242   PetscCall(MatSetLocalToGlobalMapping(j21,vISMap,eISMap));
2243   PetscCall(MatSetLocalToGlobalMapping(j22,vISMap,vISMap));
2244 
2245   PetscCall(MatSetUp(j11));
2246   PetscCall(MatSetUp(j12));
2247   PetscCall(MatSetUp(j21));
2248   PetscCall(MatSetUp(j22));
2249 
2250   PetscCall(MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J));
2251   PetscCall(MatSetUp(*J));
2252   PetscCall(MatNestSetVecType(*J,VECNEST));
2253   PetscCall(MatDestroy(&j11));
2254   PetscCall(MatDestroy(&j12));
2255   PetscCall(MatDestroy(&j21));
2256   PetscCall(MatDestroy(&j22));
2257 
2258   PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY));
2259   PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY));
2260   PetscCall(MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
2261 
2262   /* Free structures */
2263   PetscCall(ISLocalToGlobalMappingDestroy(&eISMap));
2264   PetscCall(ISLocalToGlobalMappingDestroy(&vISMap));
2265   PetscFunctionReturn(0);
2266 }
2267 
2268 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
2269 {
2270   DM_Network     *network = (DM_Network*)dm->data;
2271   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
2272   PetscInt       cstart,ncols,j,e,v;
2273   PetscBool      ghost,ghost_vc,ghost2,isNest;
2274   Mat            Juser;
2275   PetscSection   sectionGlobal;
2276   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
2277   const PetscInt *edges,*cone;
2278   MPI_Comm       comm;
2279   MatType        mtype;
2280   Vec            vd_nz,vo_nz;
2281   PetscInt       *dnnz,*onnz;
2282   PetscScalar    *vdnz,*vonz;
2283 
2284   PetscFunctionBegin;
2285   mtype = dm->mattype;
2286   PetscCall(PetscStrcmp(mtype,MATNEST,&isNest));
2287   if (isNest) {
2288     PetscCall(DMCreateMatrix_Network_Nest(dm,J));
2289     PetscCall(MatSetDM(*J,dm));
2290     PetscFunctionReturn(0);
2291   }
2292 
2293   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
2294     /* user does not provide Jacobian blocks */
2295     PetscCall(DMCreateMatrix_Plex(network->plex,J));
2296     PetscCall(MatSetDM(*J,dm));
2297     PetscFunctionReturn(0);
2298   }
2299 
2300   PetscCall(MatCreate(PetscObjectComm((PetscObject)dm),J));
2301   PetscCall(DMGetGlobalSection(network->plex,&sectionGlobal));
2302   PetscCall(PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize));
2303   PetscCall(MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE));
2304 
2305   PetscCall(MatSetType(*J,MATAIJ));
2306   PetscCall(MatSetFromOptions(*J));
2307 
2308   /* (1) Set matrix preallocation */
2309   /*------------------------------*/
2310   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
2311   PetscCall(VecCreate(comm,&vd_nz));
2312   PetscCall(VecSetSizes(vd_nz,localSize,PETSC_DECIDE));
2313   PetscCall(VecSetFromOptions(vd_nz));
2314   PetscCall(VecSet(vd_nz,0.0));
2315   PetscCall(VecDuplicate(vd_nz,&vo_nz));
2316 
2317   /* Set preallocation for edges */
2318   /*-----------------------------*/
2319   PetscCall(DMNetworkGetEdgeRange(dm,&eStart,&eEnd));
2320 
2321   PetscCall(PetscMalloc1(localSize,&rows));
2322   for (e=eStart; e<eEnd; e++) {
2323     /* Get row indices */
2324     PetscCall(DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart));
2325     PetscCall(PetscSectionGetDof(network->DofSection,e,&nrows));
2326     if (nrows) {
2327       for (j=0; j<nrows; j++) rows[j] = j + rstart;
2328 
2329       /* Set preallocation for connected vertices */
2330       PetscCall(DMNetworkGetConnectedVertices(dm,e,&cone));
2331       for (v=0; v<2; v++) {
2332         PetscCall(PetscSectionGetDof(network->DofSection,cone[v],&ncols));
2333 
2334         if (network->Je) {
2335           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
2336         } else Juser = NULL;
2337         PetscCall(DMNetworkIsGhostVertex(dm,cone[v],&ghost));
2338         PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz));
2339       }
2340 
2341       /* Set preallocation for edge self */
2342       cstart = rstart;
2343       if (network->Je) {
2344         Juser = network->Je[3*e]; /* Jacobian(e,e) */
2345       } else Juser = NULL;
2346       PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz));
2347     }
2348   }
2349 
2350   /* Set preallocation for vertices */
2351   /*--------------------------------*/
2352   PetscCall(DMNetworkGetVertexRange(dm,&vStart,&vEnd));
2353   if (vEnd - vStart) vptr = network->Jvptr;
2354 
2355   for (v=vStart; v<vEnd; v++) {
2356     /* Get row indices */
2357     PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart));
2358     PetscCall(PetscSectionGetDof(network->DofSection,v,&nrows));
2359     if (!nrows) continue;
2360 
2361     PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost));
2362     if (ghost) {
2363       PetscCall(PetscMalloc1(nrows,&rows_v));
2364     } else {
2365       rows_v = rows;
2366     }
2367 
2368     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
2369 
2370     /* Get supporting edges and connected vertices */
2371     PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges));
2372 
2373     for (e=0; e<nedges; e++) {
2374       /* Supporting edges */
2375       PetscCall(DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart));
2376       PetscCall(PetscSectionGetDof(network->DofSection,edges[e],&ncols));
2377 
2378       if (network->Jv) {
2379         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
2380       } else Juser = NULL;
2381       PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz));
2382 
2383       /* Connected vertices */
2384       PetscCall(DMNetworkGetConnectedVertices(dm,edges[e],&cone));
2385       vc = (v == cone[0]) ? cone[1]:cone[0];
2386       PetscCall(DMNetworkIsGhostVertex(dm,vc,&ghost_vc));
2387 
2388       PetscCall(PetscSectionGetDof(network->DofSection,vc,&ncols));
2389 
2390       if (network->Jv) {
2391         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
2392       } else Juser = NULL;
2393       if (ghost_vc||ghost) {
2394         ghost2 = PETSC_TRUE;
2395       } else {
2396         ghost2 = PETSC_FALSE;
2397       }
2398       PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz));
2399     }
2400 
2401     /* Set preallocation for vertex self */
2402     PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost));
2403     if (!ghost) {
2404       PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart));
2405       if (network->Jv) {
2406         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
2407       } else Juser = NULL;
2408       PetscCall(MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz));
2409     }
2410     if (ghost) {
2411       PetscCall(PetscFree(rows_v));
2412     }
2413   }
2414 
2415   PetscCall(VecAssemblyBegin(vd_nz));
2416   PetscCall(VecAssemblyBegin(vo_nz));
2417 
2418   PetscCall(PetscMalloc2(localSize,&dnnz,localSize,&onnz));
2419 
2420   PetscCall(VecAssemblyEnd(vd_nz));
2421   PetscCall(VecAssemblyEnd(vo_nz));
2422 
2423   PetscCall(VecGetArray(vd_nz,&vdnz));
2424   PetscCall(VecGetArray(vo_nz,&vonz));
2425   for (j=0; j<localSize; j++) {
2426     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
2427     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
2428   }
2429   PetscCall(VecRestoreArray(vd_nz,&vdnz));
2430   PetscCall(VecRestoreArray(vo_nz,&vonz));
2431   PetscCall(VecDestroy(&vd_nz));
2432   PetscCall(VecDestroy(&vo_nz));
2433 
2434   PetscCall(MatSeqAIJSetPreallocation(*J,0,dnnz));
2435   PetscCall(MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz));
2436   PetscCall(MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE));
2437 
2438   PetscCall(PetscFree2(dnnz,onnz));
2439 
2440   /* (2) Set matrix entries for edges */
2441   /*----------------------------------*/
2442   for (e=eStart; e<eEnd; e++) {
2443     /* Get row indices */
2444     PetscCall(DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart));
2445     PetscCall(PetscSectionGetDof(network->DofSection,e,&nrows));
2446     if (nrows) {
2447       for (j=0; j<nrows; j++) rows[j] = j + rstart;
2448 
2449       /* Set matrix entries for connected vertices */
2450       PetscCall(DMNetworkGetConnectedVertices(dm,e,&cone));
2451       for (v=0; v<2; v++) {
2452         PetscCall(DMNetworkGetGlobalVecOffset(dm,cone[v],ALL_COMPONENTS,&cstart));
2453         PetscCall(PetscSectionGetDof(network->DofSection,cone[v],&ncols));
2454 
2455         if (network->Je) {
2456           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
2457         } else Juser = NULL;
2458         PetscCall(MatSetblock_private(Juser,nrows,rows,ncols,cstart,J));
2459       }
2460 
2461       /* Set matrix entries for edge self */
2462       cstart = rstart;
2463       if (network->Je) {
2464         Juser = network->Je[3*e]; /* Jacobian(e,e) */
2465       } else Juser = NULL;
2466       PetscCall(MatSetblock_private(Juser,nrows,rows,nrows,cstart,J));
2467     }
2468   }
2469 
2470   /* Set matrix entries for vertices */
2471   /*---------------------------------*/
2472   for (v=vStart; v<vEnd; v++) {
2473     /* Get row indices */
2474     PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart));
2475     PetscCall(PetscSectionGetDof(network->DofSection,v,&nrows));
2476     if (!nrows) continue;
2477 
2478     PetscCall(DMNetworkIsGhostVertex(dm,v,&ghost));
2479     if (ghost) {
2480       PetscCall(PetscMalloc1(nrows,&rows_v));
2481     } else {
2482       rows_v = rows;
2483     }
2484     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
2485 
2486     /* Get supporting edges and connected vertices */
2487     PetscCall(DMNetworkGetSupportingEdges(dm,v,&nedges,&edges));
2488 
2489     for (e=0; e<nedges; e++) {
2490       /* Supporting edges */
2491       PetscCall(DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart));
2492       PetscCall(PetscSectionGetDof(network->DofSection,edges[e],&ncols));
2493 
2494       if (network->Jv) {
2495         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
2496       } else Juser = NULL;
2497       PetscCall(MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J));
2498 
2499       /* Connected vertices */
2500       PetscCall(DMNetworkGetConnectedVertices(dm,edges[e],&cone));
2501       vc = (v == cone[0]) ? cone[1]:cone[0];
2502 
2503       PetscCall(DMNetworkGetGlobalVecOffset(dm,vc,ALL_COMPONENTS,&cstart));
2504       PetscCall(PetscSectionGetDof(network->DofSection,vc,&ncols));
2505 
2506       if (network->Jv) {
2507         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
2508       } else Juser = NULL;
2509       PetscCall(MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J));
2510     }
2511 
2512     /* Set matrix entries for vertex self */
2513     if (!ghost) {
2514       PetscCall(DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart));
2515       if (network->Jv) {
2516         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
2517       } else Juser = NULL;
2518       PetscCall(MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J));
2519     }
2520     if (ghost) {
2521       PetscCall(PetscFree(rows_v));
2522     }
2523   }
2524   PetscCall(PetscFree(rows));
2525 
2526   PetscCall(MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY));
2527   PetscCall(MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY));
2528 
2529   PetscCall(MatSetDM(*J,dm));
2530   PetscFunctionReturn(0);
2531 }
2532 
2533 PetscErrorCode DMDestroy_Network(DM dm)
2534 {
2535   DM_Network     *network = (DM_Network*)dm->data;
2536   PetscInt       j,np;
2537 
2538   PetscFunctionBegin;
2539   if (--network->refct > 0) PetscFunctionReturn(0);
2540   PetscCall(PetscFree(network->Je));
2541   if (network->Jv) {
2542     PetscCall(PetscFree(network->Jvptr));
2543     PetscCall(PetscFree(network->Jv));
2544   }
2545 
2546   PetscCall(ISLocalToGlobalMappingDestroy(&network->vertex.mapping));
2547   PetscCall(PetscSectionDestroy(&network->vertex.DofSection));
2548   PetscCall(PetscSectionDestroy(&network->vertex.GlobalDofSection));
2549   PetscCall(PetscFree(network->vltog));
2550   PetscCall(PetscSFDestroy(&network->vertex.sf));
2551   /* edge */
2552   PetscCall(ISLocalToGlobalMappingDestroy(&network->edge.mapping));
2553   PetscCall(PetscSectionDestroy(&network->edge.DofSection));
2554   PetscCall(PetscSectionDestroy(&network->edge.GlobalDofSection));
2555   PetscCall(PetscSFDestroy(&network->edge.sf));
2556   PetscCall(DMDestroy(&network->plex));
2557   PetscCall(PetscSectionDestroy(&network->DataSection));
2558   PetscCall(PetscSectionDestroy(&network->DofSection));
2559 
2560   for (j=0; j<network->Nsvtx; j++) PetscCall(PetscFree(network->svtx[j].sv));
2561   PetscCall(PetscFree(network->svtx));
2562   PetscCall(PetscFree2(network->subnetedge,network->subnetvtx));
2563 
2564   PetscCall(PetscTableDestroy(&network->svtable));
2565   PetscCall(PetscFree(network->subnet));
2566   PetscCall(PetscFree(network->component));
2567   PetscCall(PetscFree(network->componentdataarray));
2568 
2569   if (network->header) {
2570     np = network->pEnd - network->pStart;
2571     for (j=0; j < np; j++) {
2572       PetscCall(PetscFree5(network->header[j].size,network->header[j].key,network->header[j].offset,network->header[j].nvar,network->header[j].offsetvarrel));
2573       PetscCall(PetscFree(network->cvalue[j].data));
2574     }
2575     PetscCall(PetscFree2(network->header,network->cvalue));
2576   }
2577   PetscCall(PetscFree(network));
2578   PetscFunctionReturn(0);
2579 }
2580 
2581 PetscErrorCode DMView_Network(DM dm,PetscViewer viewer)
2582 {
2583   PetscBool      iascii;
2584   PetscMPIInt    rank;
2585 
2586   PetscFunctionBegin;
2587   PetscValidHeaderSpecific(dm,DM_CLASSID, 1);
2588   PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2);
2589   PetscCallMPI(MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank));
2590   PetscCall(PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii));
2591   if (iascii) {
2592     const PetscInt *cone,*vtx,*edges;
2593     PetscInt       vfrom,vto,i,j,nv,ne,nsv,p,nsubnet;
2594     DM_Network     *network = (DM_Network*)dm->data;
2595 
2596     nsubnet = network->Nsubnet; /* num of subnetworks */
2597     if (rank == 0) {
2598       PetscCall(PetscPrintf(PETSC_COMM_SELF,"  NSubnets: %" PetscInt_FMT "; NEdges: %" PetscInt_FMT "; NVertices: %" PetscInt_FMT "; NSharedVertices: %" PetscInt_FMT ".\n",nsubnet,network->NEdges,network->NVertices,network->Nsvtx));
2599     }
2600 
2601     PetscCall(DMNetworkGetSharedVertices(dm,&nsv,NULL));
2602     PetscCall(PetscViewerASCIIPushSynchronized(viewer));
2603     PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "  [%d] nEdges: %" PetscInt_FMT "; nVertices: %" PetscInt_FMT "; nSharedVertices: %" PetscInt_FMT "\n",rank,network->nEdges,network->nVertices,nsv));
2604 
2605     for (i=0; i<nsubnet; i++) {
2606       PetscCall(DMNetworkGetSubnetwork(dm,i,&nv,&ne,&vtx,&edges));
2607       if (ne) {
2608         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "     Subnet %" PetscInt_FMT ": nEdges %" PetscInt_FMT ", nVertices(include shared vertices) %" PetscInt_FMT "\n",i,ne,nv));
2609         for (j=0; j<ne; j++) {
2610           p = edges[j];
2611           PetscCall(DMNetworkGetConnectedVertices(dm,p,&cone));
2612           PetscCall(DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom));
2613           PetscCall(DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto));
2614           PetscCall(DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p));
2615           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "       edge %" PetscInt_FMT ": %" PetscInt_FMT " ----> %" PetscInt_FMT "\n",p,vfrom,vto));
2616         }
2617       }
2618     }
2619 
2620     /* Shared vertices */
2621     PetscCall(DMNetworkGetSharedVertices(dm,NULL,&vtx));
2622     if (nsv) {
2623       PetscInt       gidx;
2624       PetscBool      ghost;
2625       const PetscInt *sv=NULL;
2626 
2627       PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "     SharedVertices:\n"));
2628       for (i=0; i<nsv; i++) {
2629         PetscCall(DMNetworkIsGhostVertex(dm,vtx[i],&ghost));
2630         if (ghost) continue;
2631 
2632         PetscCall(DMNetworkSharedVertexGetInfo(dm,vtx[i],&gidx,&nv,&sv));
2633         PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "       svtx %" PetscInt_FMT ": global index %" PetscInt_FMT ", subnet[%" PetscInt_FMT "].%" PetscInt_FMT " ---->\n",i,gidx,sv[0],sv[1]));
2634         for (j=1; j<nv; j++) {
2635           PetscCall(PetscViewerASCIISynchronizedPrintf(viewer, "                                           ----> subnet[%" PetscInt_FMT "].%" PetscInt_FMT "\n",sv[2*j],sv[2*j+1]));
2636         }
2637       }
2638     }
2639     PetscCall(PetscViewerFlush(viewer));
2640     PetscCall(PetscViewerASCIIPopSynchronized(viewer));
2641   } else PetscCheck(iascii,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"Viewer type %s not yet supported for DMNetwork writing",((PetscObject)viewer)->type_name);
2642   PetscFunctionReturn(0);
2643 }
2644 
2645 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
2646 {
2647   DM_Network     *network = (DM_Network*)dm->data;
2648 
2649   PetscFunctionBegin;
2650   PetscCall(DMGlobalToLocalBegin(network->plex,g,mode,l));
2651   PetscFunctionReturn(0);
2652 }
2653 
2654 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
2655 {
2656   DM_Network     *network = (DM_Network*)dm->data;
2657 
2658   PetscFunctionBegin;
2659   PetscCall(DMGlobalToLocalEnd(network->plex,g,mode,l));
2660   PetscFunctionReturn(0);
2661 }
2662 
2663 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
2664 {
2665   DM_Network     *network = (DM_Network*)dm->data;
2666 
2667   PetscFunctionBegin;
2668   PetscCall(DMLocalToGlobalBegin(network->plex,l,mode,g));
2669   PetscFunctionReturn(0);
2670 }
2671 
2672 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
2673 {
2674   DM_Network     *network = (DM_Network*)dm->data;
2675 
2676   PetscFunctionBegin;
2677   PetscCall(DMLocalToGlobalEnd(network->plex,l,mode,g));
2678   PetscFunctionReturn(0);
2679 }
2680 
2681 /*@
2682   DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index
2683 
2684   Not collective
2685 
2686   Input Parameters:
2687 + dm - the dm object
2688 - vloc - local vertex ordering, start from 0
2689 
2690   Output Parameters:
2691 .  vg  - global vertex ordering, start from 0
2692 
2693   Level: advanced
2694 
2695 .seealso: DMNetworkSetVertexLocalToGlobalOrdering()
2696 @*/
2697 PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm,PetscInt vloc,PetscInt *vg)
2698 {
2699   DM_Network  *network = (DM_Network*)dm->data;
2700   PetscInt    *vltog = network->vltog;
2701 
2702   PetscFunctionBegin;
2703   PetscCheck(vltog,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkSetVertexLocalToGlobalOrdering() first");
2704   *vg = vltog[vloc];
2705   PetscFunctionReturn(0);
2706 }
2707 
2708 /*@
2709   DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map
2710 
2711   Collective
2712 
2713   Input Parameters:
2714 . dm - the dm object
2715 
2716   Level: advanced
2717 
2718 .seealso: DMNetworkGetGlobalVertexIndex()
2719 @*/
2720 PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm)
2721 {
2722   DM_Network        *network=(DM_Network*)dm->data;
2723   MPI_Comm          comm;
2724   PetscMPIInt       rank,size,*displs=NULL,*recvcounts=NULL,remoterank;
2725   PetscBool         ghost;
2726   PetscInt          *vltog,nroots,nleaves,i,*vrange,k,N,lidx;
2727   const PetscSFNode *iremote;
2728   PetscSF           vsf;
2729   Vec               Vleaves,Vleaves_seq;
2730   VecScatter        ctx;
2731   PetscScalar       *varr,val;
2732   const PetscScalar *varr_read;
2733 
2734   PetscFunctionBegin;
2735   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
2736   PetscCallMPI(MPI_Comm_size(comm,&size));
2737   PetscCallMPI(MPI_Comm_rank(comm,&rank));
2738 
2739   if (size == 1) {
2740     nroots = network->vEnd - network->vStart;
2741     PetscCall(PetscMalloc1(nroots, &vltog));
2742     for (i=0; i<nroots; i++) vltog[i] = i;
2743     network->vltog = vltog;
2744     PetscFunctionReturn(0);
2745   }
2746 
2747   PetscCheck(network->distributecalled,comm, PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkDistribute() first");
2748   if (network->vltog) {
2749     PetscCall(PetscFree(network->vltog));
2750   }
2751 
2752   PetscCall(DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping));
2753   PetscCall(PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf));
2754   vsf = network->vertex.sf;
2755 
2756   PetscCall(PetscMalloc3(size+1,&vrange,size,&displs,size,&recvcounts));
2757   PetscCall(PetscSFGetGraph(vsf,&nroots,&nleaves,NULL,&iremote));
2758 
2759   for (i=0; i<size; i++) { displs[i] = i; recvcounts[i] = 1;}
2760 
2761   i         = nroots - nleaves; /* local number of vertices, excluding ghosts */
2762   vrange[0] = 0;
2763   PetscCallMPI(MPI_Allgatherv(&i,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm));
2764   for (i=2; i<size+1; i++) {vrange[i] += vrange[i-1];}
2765 
2766   PetscCall(PetscMalloc1(nroots, &vltog));
2767   network->vltog = vltog;
2768 
2769   /* Set vltog for non-ghost vertices */
2770   k = 0;
2771   for (i=0; i<nroots; i++) {
2772     PetscCall(DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost));
2773     if (ghost) continue;
2774     vltog[i] = vrange[rank] + k++;
2775   }
2776   PetscCall(PetscFree3(vrange,displs,recvcounts));
2777 
2778   /* Set vltog for ghost vertices */
2779   /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */
2780   PetscCall(VecCreate(comm,&Vleaves));
2781   PetscCall(VecSetSizes(Vleaves,2*nleaves,PETSC_DETERMINE));
2782   PetscCall(VecSetFromOptions(Vleaves));
2783   PetscCall(VecGetArray(Vleaves,&varr));
2784   for (i=0; i<nleaves; i++) {
2785     varr[2*i]   = (PetscScalar)(iremote[i].rank);  /* rank of remote process */
2786     varr[2*i+1] = (PetscScalar)(iremote[i].index); /* local index in remote process */
2787   }
2788   PetscCall(VecRestoreArray(Vleaves,&varr));
2789 
2790   /* (b) scatter local info to remote processes via VecScatter() */
2791   PetscCall(VecScatterCreateToAll(Vleaves,&ctx,&Vleaves_seq));
2792   PetscCall(VecScatterBegin(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD));
2793   PetscCall(VecScatterEnd(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD));
2794 
2795   /* (c) convert local indices to global indices in parallel vector Vleaves */
2796   PetscCall(VecGetSize(Vleaves_seq,&N));
2797   PetscCall(VecGetArrayRead(Vleaves_seq,&varr_read));
2798   for (i=0; i<N; i+=2) {
2799     remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]);
2800     if (remoterank == rank) {
2801       k = i+1; /* row number */
2802       lidx = (PetscInt)PetscRealPart(varr_read[i+1]);
2803       val  = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */
2804       PetscCall(VecSetValues(Vleaves,1,&k,&val,INSERT_VALUES));
2805     }
2806   }
2807   PetscCall(VecRestoreArrayRead(Vleaves_seq,&varr_read));
2808   PetscCall(VecAssemblyBegin(Vleaves));
2809   PetscCall(VecAssemblyEnd(Vleaves));
2810 
2811   /* (d) Set vltog for ghost vertices by copying local values of Vleaves */
2812   PetscCall(VecGetArrayRead(Vleaves,&varr_read));
2813   k = 0;
2814   for (i=0; i<nroots; i++) {
2815     PetscCall(DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost));
2816     if (!ghost) continue;
2817     vltog[i] = (PetscInt)PetscRealPart(varr_read[2*k+1]); k++;
2818   }
2819   PetscCall(VecRestoreArrayRead(Vleaves,&varr_read));
2820 
2821   PetscCall(VecDestroy(&Vleaves));
2822   PetscCall(VecDestroy(&Vleaves_seq));
2823   PetscCall(VecScatterDestroy(&ctx));
2824   PetscFunctionReturn(0);
2825 }
2826 
2827 static inline PetscErrorCode DMISAddSize_private(DM_Network *network,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *nidx)
2828 {
2829   PetscInt                 i,j,ncomps,nvar,key,offset=0;
2830   DMNetworkComponentHeader header;
2831 
2832   PetscFunctionBegin;
2833   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
2834   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
2835   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
2836 
2837   for (i=0; i<ncomps; i++) {
2838     key  = header->key[i];
2839     nvar = header->nvar[i];
2840     for (j=0; j<numkeys; j++) {
2841       if (key == keys[j]) {
2842         if (!blocksize || blocksize[j] == -1) {
2843           *nidx += nvar;
2844         } else {
2845           *nidx += nselectedvar[j]*nvar/blocksize[j];
2846         }
2847       }
2848     }
2849   }
2850   PetscFunctionReturn(0);
2851 }
2852 
2853 static inline PetscErrorCode DMISComputeIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx)
2854 {
2855   PetscInt                 i,j,ncomps,nvar,key,offsetg,k,k1,offset=0;
2856   DM_Network               *network = (DM_Network*)dm->data;
2857   DMNetworkComponentHeader header;
2858 
2859   PetscFunctionBegin;
2860   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
2861   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
2862   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
2863 
2864   for (i=0; i<ncomps; i++) {
2865     key  = header->key[i];
2866     nvar = header->nvar[i];
2867     for (j=0; j<numkeys; j++) {
2868       if (key != keys[j]) continue;
2869 
2870       PetscCall(DMNetworkGetGlobalVecOffset(dm,p,i,&offsetg));
2871       if (!blocksize || blocksize[j] == -1) {
2872         for (k=0; k<nvar; k++) idx[(*ii)++] = offsetg + k;
2873       } else {
2874         for (k=0; k<nvar; k+=blocksize[j]) {
2875           for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1];
2876         }
2877       }
2878     }
2879   }
2880   PetscFunctionReturn(0);
2881 }
2882 
2883 /*@
2884   DMNetworkCreateIS - Create an index set object from the global vector of the network
2885 
2886   Collective
2887 
2888   Input Parameters:
2889 + dm - DMNetwork object
2890 . numkeys - number of keys
2891 . keys - array of keys that define the components of the variables you wish to extract
2892 . blocksize - block size of the variables associated to the component
2893 . nselectedvar - number of variables in each block to select
2894 - selectedvar - the offset into the block of each variable in each block to select
2895 
2896   Output Parameters:
2897 . is - the index set
2898 
2899   Level: Advanced
2900 
2901   Notes:
2902     Use blocksize[i] of -1 to indicate select all the variables of the i-th component, for which nselectvar[i] and selectedvar[i] are ignored. Use NULL, NULL, NULL to indicate for all selected components one wishes to obtain all the values of that component. For example, DMNetworkCreateIS(dm,1,&key,NULL,NULL,NULL,&is) will return an is that extracts all the variables for the 0-th component.
2903 
2904 .seealso: DMNetworkCreate(), ISCreateGeneral(), DMNetworkCreateLocalIS()
2905 @*/
2906 PetscErrorCode DMNetworkCreateIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is)
2907 {
2908   MPI_Comm       comm;
2909   DM_Network     *network = (DM_Network*)dm->data;
2910   PetscInt       i,p,estart,eend,vstart,vend,nidx,*idx;
2911   PetscBool      ghost;
2912 
2913   PetscFunctionBegin;
2914   PetscCall(PetscObjectGetComm((PetscObject)dm,&comm));
2915 
2916   /* Check input parameters */
2917   for (i=0; i<numkeys; i++) {
2918     if (!blocksize || blocksize[i] == -1) continue;
2919     PetscCheck(nselectedvar[i] <= blocksize[i],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"number of selectedvariables %" PetscInt_FMT " cannot be larger than blocksize %" PetscInt_FMT "",nselectedvar[i],blocksize[i]);
2920   }
2921 
2922   PetscCall(DMNetworkGetEdgeRange(dm,&estart,&eend));
2923   PetscCall(DMNetworkGetVertexRange(dm,&vstart,&vend));
2924 
2925   /* Get local number of idx */
2926   nidx = 0;
2927   for (p=estart; p<eend; p++) {
2928     PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx));
2929   }
2930   for (p=vstart; p<vend; p++) {
2931     PetscCall(DMNetworkIsGhostVertex(dm,p,&ghost));
2932     if (ghost) continue;
2933     PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx));
2934   }
2935 
2936   /* Compute idx */
2937   PetscCall(PetscMalloc1(nidx,&idx));
2938   i = 0;
2939   for (p=estart; p<eend; p++) {
2940     PetscCall(DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx));
2941   }
2942   for (p=vstart; p<vend; p++) {
2943     PetscCall(DMNetworkIsGhostVertex(dm,p,&ghost));
2944     if (ghost) continue;
2945     PetscCall(DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx));
2946   }
2947 
2948   /* Create is */
2949   PetscCall(ISCreateGeneral(comm,nidx,idx,PETSC_COPY_VALUES,is));
2950   PetscCall(PetscFree(idx));
2951   PetscFunctionReturn(0);
2952 }
2953 
2954 static inline PetscErrorCode DMISComputeLocalIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx)
2955 {
2956   PetscInt                 i,j,ncomps,nvar,key,offsetl,k,k1,offset=0;
2957   DM_Network               *network = (DM_Network*)dm->data;
2958   DMNetworkComponentHeader header;
2959 
2960   PetscFunctionBegin;
2961   PetscCall(PetscSectionGetOffset(network->DataSection,p,&offset));
2962   ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
2963   header = (DMNetworkComponentHeader)(network->componentdataarray+offset);
2964 
2965   for (i=0; i<ncomps; i++) {
2966     key  = header->key[i];
2967     nvar = header->nvar[i];
2968     for (j=0; j<numkeys; j++) {
2969       if (key != keys[j]) continue;
2970 
2971       PetscCall(DMNetworkGetLocalVecOffset(dm,p,i,&offsetl));
2972       if (!blocksize || blocksize[j] == -1) {
2973         for (k=0; k<nvar; k++) idx[(*ii)++] = offsetl + k;
2974       } else {
2975         for (k=0; k<nvar; k+=blocksize[j]) {
2976           for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1];
2977         }
2978       }
2979     }
2980   }
2981   PetscFunctionReturn(0);
2982 }
2983 
2984 /*@
2985   DMNetworkCreateLocalIS - Create an index set object from the local vector of the network
2986 
2987   Not collective
2988 
2989   Input Parameters:
2990 + dm - DMNetwork object
2991 . numkeys - number of keys
2992 . keys - array of keys that define the components of the variables you wish to extract
2993 . blocksize - block size of the variables associated to the component
2994 . nselectedvar - number of variables in each block to select
2995 - selectedvar - the offset into the block of each variable in each block to select
2996 
2997   Output Parameters:
2998 . is - the index set
2999 
3000   Level: Advanced
3001 
3002   Notes:
3003     Use blocksize[i] of -1 to indicate select all the variables of the i-th component, for which nselectvar[i] and selectedvar[i] are ignored. Use NULL, NULL, NULL to indicate for all selected components one wishes to obtain all the values of that component. For example, DMNetworkCreateLocalIS(dm,1,&key,NULL,NULL,NULL,&is) will return an is that extracts all the variables for the 0-th component.
3004 
3005 .seealso: DMNetworkCreate(), DMNetworkCreateIS, ISCreateGeneral()
3006 @*/
3007 PetscErrorCode DMNetworkCreateLocalIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is)
3008 {
3009   DM_Network     *network = (DM_Network*)dm->data;
3010   PetscInt       i,p,pstart,pend,nidx,*idx;
3011 
3012   PetscFunctionBegin;
3013   /* Check input parameters */
3014   for (i=0; i<numkeys; i++) {
3015     if (!blocksize || blocksize[i] == -1) continue;
3016     PetscCheck(nselectedvar[i] <= blocksize[i],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"number of selectedvariables %" PetscInt_FMT " cannot be larger than blocksize %" PetscInt_FMT "",nselectedvar[i],blocksize[i]);
3017   }
3018 
3019   pstart = network->pStart;
3020   pend   = network->pEnd;
3021 
3022   /* Get local number of idx */
3023   nidx = 0;
3024   for (p=pstart; p<pend; p++) {
3025     PetscCall(DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx));
3026   }
3027 
3028   /* Compute local idx */
3029   PetscCall(PetscMalloc1(nidx,&idx));
3030   i = 0;
3031   for (p=pstart; p<pend; p++) {
3032     PetscCall(DMISComputeLocalIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx));
3033   }
3034 
3035   /* Create is */
3036   PetscCall(ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,is));
3037   PetscCall(PetscFree(idx));
3038   PetscFunctionReturn(0);
3039 }
3040