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