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