xref: /petsc/src/dm/impls/network/network.c (revision 7765340c1d68d7d22cd343b2706986f68198c29a)
1 #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
2 #include <petscdmplex.h>
3 #include <petscsf.h>
4 
5 /*@
6   DMNetworkGetPlex - Gets the Plex DM associated with this network DM
7 
8   Not collective
9 
10   Input Parameters:
11 + netdm - the dm object
12 - plexmdm - the plex dm object
13 
14   Level: Advanced
15 
16 .seealso: DMNetworkCreate()
17 @*/
18 PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm)
19 {
20   DM_Network     *network = (DM_Network*) netdm->data;
21 
22   PetscFunctionBegin;
23   *plexdm = network->plex;
24   PetscFunctionReturn(0);
25 }
26 
27 /*@
28   DMNetworkSetSizes - Sets the number of subnetworks,local and global vertices and edges for each subnetwork.
29 
30   Collective on DM
31 
32   Input Parameters:
33 + dm - the dm object
34 . Nsubnet - number of subnetworks
35 . nV - number of local vertices for each subnetwork
36 . nE - number of local edges for each subnetwork
37 . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork
38 - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork
39 
40    Notes
41    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
42 
43    You cannot change the sizes once they have been set
44 
45    Level: intermediate
46 
47 .seealso: DMNetworkCreate()
48 @*/
49 PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt Nsubnet, PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[])
50 {
51   PetscErrorCode ierr;
52   DM_Network     *network = (DM_Network*) dm->data;
53   PetscInt       a[2],b[2],i;
54 
55   PetscFunctionBegin;
56   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
57   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
58 
59   if(Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2);
60   if(network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
61 
62   network->nsubnet = Nsubnet;
63   ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr);
64   for(i=0; i < network->nsubnet; i++) {
65     if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],5);
66     if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],6);
67     if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]);
68     if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]);
69     a[0] = nV[i]; a[1] = nE[i];
70     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
71     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];
72 
73     network->subnet[i].id = i;
74 
75     network->subnet[i].nvtx = nV[i];
76     network->subnet[i].vStart = network->nVertices;
77     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
78     network->nVertices += network->subnet[i].nvtx;
79     network->NVertices += network->subnet[i].Nvtx;
80 
81     network->subnet[i].nedge = nE[i];
82     network->subnet[i].eStart = network->nEdges;
83     network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge;
84     network->nEdges += network->subnet[i].nedge;
85     network->NEdges += network->subnet[i].Nedge;
86   }
87   PetscFunctionReturn(0);
88 }
89 
90 PetscErrorCode DMNetworkSetSizesCoupled(DM dm,PetscInt Nsubnet,PetscInt NsubnetCouple,PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[])
91 {
92   PetscErrorCode ierr;
93   DM_Network     *network = (DM_Network*) dm->data;
94   PetscInt       a[2],b[2],i;
95 
96   PetscFunctionBegin;
97   printf("DMNetworkSetSizesCoupled...Nsubnet %d, NsubnetCouple %d\n",Nsubnet,NsubnetCouple);
98   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
99   if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet);
100   if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple);
101 
102   if (Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2);
103   if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network");
104 
105   network->nsubnet  = Nsubnet + NsubnetCouple;
106   Nsubnet          += NsubnetCouple;
107   network->ncsubnet = NsubnetCouple;
108   ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr);
109   for(i=0; i < Nsubnet; i++) {
110     if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],6);
111     if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],7);
112     if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]);
113     if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]);
114     a[0] = nV[i]; a[1] = nE[i];
115     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
116     network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1];
117 
118     network->subnet[i].id = i;
119 
120     network->subnet[i].nvtx = nV[i];
121     network->subnet[i].vStart = network->nVertices;
122     network->subnet[i].vEnd   = network->subnet[i].vStart + network->subnet[i].Nvtx;
123     network->nVertices += network->subnet[i].nvtx;
124     network->NVertices += network->subnet[i].Nvtx;
125 
126     network->subnet[i].nedge = nE[i];
127     network->subnet[i].eStart = network->nEdges;
128     network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge;
129     network->nEdges += network->subnet[i].nedge;
130     network->NEdges += network->subnet[i].Nedge;
131   }
132   PetscFunctionReturn(0);
133 }
134 
135 /*@
136   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
137 
138   Logically collective on DM
139 
140   Input Parameters:
141 . edges - list of edges for each subnetwork
142 
143   Notes:
144   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
145   not be destroyed before the call to DMNetworkLayoutSetUp
146 
147   Level: intermediate
148 
149 .seealso: DMNetworkCreate, DMNetworkSetSizes
150 @*/
151 PetscErrorCode DMNetworkSetEdgeList(DM dm, int *edgelist[])
152 {
153   DM_Network *network = (DM_Network*) dm->data;
154   PetscInt   i;
155 
156   PetscFunctionBegin;
157   for(i=0; i < network->nsubnet; i++) {
158     network->subnet[i].edgelist = edgelist[i];
159   }
160   PetscFunctionReturn(0);
161 }
162 
163 /*@
164   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
165 
166   Collective on DM
167 
168   Input Parameters
169 . DM - the dmnetwork object
170 
171   Notes:
172   This routine should be called after the network sizes and edgelists have been provided. It creates
173   the bare layout of the network and sets up the network to begin insertion of components.
174 
175   All the components should be registered before calling this routine.
176 
177   Level: intermediate
178 
179 .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
180 @*/
181 PetscErrorCode DMNetworkLayoutSetUp(DM dm)
182 {
183   PetscErrorCode ierr;
184   DM_Network     *network = (DM_Network*) dm->data;
185   PetscInt       dim = 1; /* One dimensional network */
186   PetscInt       numCorners=2;
187   PetscInt       spacedim=2;
188   double         *vertexcoords=NULL;
189   PetscInt       i,j;
190   PetscInt       ndata;
191   PetscInt       ctr=0;
192 
193   PetscFunctionBegin;
194   if (network->nVertices) {
195     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
196   }
197 
198   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
199   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
200   for(i=0; i < network->nsubnet; i++) {
201     for(j = 0; j < network->subnet[i].nedge; j++) {
202       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
203       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
204       ctr++;
205     }
206   }
207 
208 #if 0
209   for(i=0; i < network->nEdges; i++) {
210     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
211   }
212 #endif
213 
214   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
215   if (network->nVertices) {
216     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
217   }
218   ierr = PetscFree(network->edges);CHKERRQ(ierr);
219 
220   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
221   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
222   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
223 
224   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
225   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
226   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
227   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
228 
229   /* Create vertices and edges array for the subnetworks */
230   for(j=0; j < network->nsubnet; j++) {
231     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
232     ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr);
233     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
234        These get updated when the vertices and edges are added. */
235     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
236   }
237 
238   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
239   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
240   for(i=network->eStart; i < network->eEnd; i++) {
241     network->header[i].index = i;   /* Global edge number */
242     for(j=0; j < network->nsubnet; j++) {
243       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
244 	network->header[i].subnetid = j; /* Subnetwork id */
245 	network->subnet[j].edges[network->subnet[j].nedge++] = i;
246 	break;
247       }
248     }
249     network->header[i].ndata = 0;
250     ndata = network->header[i].ndata;
251     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
252     network->header[i].offset[ndata] = 0;
253   }
254 
255   for(i=network->vStart; i < network->vEnd; i++) {
256     network->header[i].index = i - network->vStart;
257     for(j=0; j < network->nsubnet; j++) {
258       if((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
259 	network->header[i].subnetid = j;
260 	network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
261 	break;
262       }
263     }
264     network->header[i].ndata = 0;
265     ndata = network->header[i].ndata;
266     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
267     network->header[i].offset[ndata] = 0;
268   }
269 
270   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
271   PetscFunctionReturn(0);
272 }
273 
274 PetscErrorCode DMNetworkLayoutSetUpCoupled(DM dm)
275 {
276   PetscErrorCode ierr;
277   DM_Network     *network = (DM_Network*)dm->data;
278   PetscInt       dim = 1; /* One dimensional network */
279   PetscInt       numCorners=2,spacedim=2;
280   double         *vertexcoords=NULL;
281   PetscInt       i,j,ndata,ctr=0,nsubnet;
282   PetscInt       *edgelist_couple=NULL,k,netid,vid;
283 
284   PetscFunctionBegin;
285   if (network->nVertices) {
286     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
287   }
288 
289   /* Create the edgelist for the network by concatenating edgelists of the subnetworks */
290   nsubnet = network->nsubnet - network->ncsubnet;
291   ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr);
292   for (i=0; i < nsubnet; i++) {
293     for (j = 0; j < network->subnet[i].nedge; j++) {
294       network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j];
295       network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1];
296       ctr++;
297     }
298   }
299 
300   i       = nsubnet; /* netid of coupling subnet */
301   nsubnet = network->nsubnet;
302   while (i < nsubnet) {
303     edgelist_couple = network->subnet[i].edgelist;
304     k = 0;
305     for (j = 0; j < network->subnet[i].nedge; j++) {
306       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
307       network->edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2;
308 
309       netid = edgelist_couple[k]; vid = edgelist_couple[k+1];
310       network->edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2;
311       ctr++;
312     }
313     i++;
314   }
315 
316 #if 1
317   for(i=0; i < network->nEdges; i++) {
318     ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr);
319     printf("\n");
320   }
321 #endif
322 
323   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
324   if (network->nVertices) {
325     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
326   }
327   ierr = PetscFree(network->edges);CHKERRQ(ierr);
328 
329   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
330   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
331   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
332 
333   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
334   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
335   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
336   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
337 
338   /* Create vertices and edges array for the subnetworks */
339   for (j=0; j < network->nsubnet; j++) {
340     ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr);
341     ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr);
342     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
343        These get updated when the vertices and edges are added. */
344     network->subnet[j].nvtx = network->subnet[j].nedge = 0;
345   }
346 
347   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
348   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
349   for (i=network->eStart; i < network->eEnd; i++) {
350     network->header[i].index = i;   /* Global edge number */
351     for (j=0; j < network->nsubnet; j++) {
352       if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) {
353 	network->header[i].subnetid = j; /* Subnetwork id */
354 	network->subnet[j].edges[network->subnet[j].nedge++] = i;
355 	break;
356       }
357     }
358     network->header[i].ndata = 0;
359     ndata = network->header[i].ndata;
360     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
361     network->header[i].offset[ndata] = 0;
362   }
363 
364   for(i=network->vStart; i < network->vEnd; i++) {
365     network->header[i].index = i - network->vStart;
366     for (j=0; j < network->nsubnet; j++) {
367       if ((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) {
368 	network->header[i].subnetid = j;
369 	network->subnet[j].vertices[network->subnet[j].nvtx++] = i;
370 	break;
371       }
372     }
373     network->header[i].ndata = 0;
374     ndata = network->header[i].ndata;
375     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
376     network->header[i].offset[ndata] = 0;
377   }
378 
379   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
380   PetscFunctionReturn(0);
381 }
382 
383 /*@C
384   DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork
385 
386   Input Parameters
387 + dm   - the number object
388 - id   - the ID (integer) of the subnetwork
389 
390   Output Parameters
391 + nv    - number of vertices (local)
392 . ne    - number of edges (local)
393 . vtx   - local vertices for this subnetwork
394 . edge  - local edges for this subnetwork
395 
396   Notes:
397   Cannot call this routine before DMNetworkLayoutSetup()
398 
399 .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
400 @*/
401 PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge)
402 {
403   DM_Network     *network = (DM_Network*) netdm->data;
404 
405   PetscFunctionBegin;
406   *nv = network->subnet[id].nvtx;
407   *ne = network->subnet[id].nedge;
408   *vtx = network->subnet[id].vertices;
409   *edge = network->subnet[id].edges;
410   PetscFunctionReturn(0);
411 }
412 
413 /*@C
414   DMNetworkRegisterComponent - Registers the network component
415 
416   Logically collective on DM
417 
418   Input Parameters
419 + dm   - the network object
420 . name - the component name
421 - size - the storage size in bytes for this component data
422 
423    Output Parameters
424 .   key - an integer key that defines the component
425 
426    Notes
427    This routine should be called by all processors before calling DMNetworkLayoutSetup().
428 
429    Level: intermediate
430 
431 .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
432 @*/
433 PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
434 {
435   PetscErrorCode        ierr;
436   DM_Network            *network = (DM_Network*) dm->data;
437   DMNetworkComponent    *component=&network->component[network->ncomponent];
438   PetscBool             flg=PETSC_FALSE;
439   PetscInt              i;
440 
441   PetscFunctionBegin;
442   for (i=0; i < network->ncomponent; i++) {
443     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
444     if (flg) {
445       *key = i;
446       PetscFunctionReturn(0);
447     }
448   }
449   if(network->ncomponent == MAX_COMPONENTS) {
450     SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS);
451   }
452 
453   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
454   component->size = size/sizeof(DMNetworkComponentGenericDataType);
455   *key = network->ncomponent;
456   network->ncomponent++;
457   PetscFunctionReturn(0);
458 }
459 
460 /*@
461   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
462 
463   Not Collective
464 
465   Input Parameters:
466 + dm - The DMNetwork object
467 
468   Output Paramters:
469 + vStart - The first vertex point
470 - vEnd   - One beyond the last vertex point
471 
472   Level: intermediate
473 
474 .seealso: DMNetworkGetEdgeRange
475 @*/
476 PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
477 {
478   DM_Network     *network = (DM_Network*)dm->data;
479 
480   PetscFunctionBegin;
481   if (vStart) *vStart = network->vStart;
482   if (vEnd) *vEnd = network->vEnd;
483   PetscFunctionReturn(0);
484 }
485 
486 /*@
487   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
488 
489   Not Collective
490 
491   Input Parameters:
492 + dm - The DMNetwork object
493 
494   Output Paramters:
495 + eStart - The first edge point
496 - eEnd   - One beyond the last edge point
497 
498   Level: intermediate
499 
500 .seealso: DMNetworkGetVertexRange
501 @*/
502 PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
503 {
504   DM_Network     *network = (DM_Network*)dm->data;
505 
506   PetscFunctionBegin;
507   if (eStart) *eStart = network->eStart;
508   if (eEnd) *eEnd = network->eEnd;
509   PetscFunctionReturn(0);
510 }
511 
512 /*@
513   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
514 
515   Not Collective
516 
517   Input Parameters:
518 + dm - DMNetwork object
519 - p  - edge point
520 
521   Output Paramters:
522 . index - user global numbering for the edge
523 
524   Level: intermediate
525 
526 .seealso: DMNetworkGetGlobalVertexIndex
527 @*/
528 PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
529 {
530   PetscErrorCode    ierr;
531   DM_Network        *network = (DM_Network*)dm->data;
532   PetscInt          offsetp;
533   DMNetworkComponentHeader header;
534 
535   PetscFunctionBegin;
536   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
537   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
538   *index = header->index;
539   PetscFunctionReturn(0);
540 }
541 
542 /*@
543   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
544 
545   Not Collective
546 
547   Input Parameters:
548 + dm - DMNetwork object
549 - p  - vertex point
550 
551   Output Paramters:
552 . index - user global numbering for the vertex
553 
554   Level: intermediate
555 
556 .seealso: DMNetworkGetGlobalEdgeIndex
557 @*/
558 PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
559 {
560   PetscErrorCode    ierr;
561   DM_Network        *network = (DM_Network*)dm->data;
562   PetscInt          offsetp;
563   DMNetworkComponentHeader header;
564 
565   PetscFunctionBegin;
566   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
567   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
568   *index = header->index;
569   PetscFunctionReturn(0);
570 }
571 
572 /*
573   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
574                                     component value from the component data array
575 
576   Not Collective
577 
578   Input Parameters:
579 + dm      - The DMNetwork object
580 . p       - vertex/edge point
581 - compnum - component number
582 
583   Output Parameters:
584 + compkey - the key obtained when registering the component
585 - offset  - offset into the component data array associated with the vertex/edge point
586 
587   Notes:
588   Typical usage:
589 
590   DMNetworkGetComponentDataArray(dm, &arr);
591   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
592   Loop over vertices or edges
593     DMNetworkGetNumComponents(dm,v,&numcomps);
594     Loop over numcomps
595       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
596       compdata = (UserCompDataType)(arr+offset);
597 
598   Level: intermediate
599 
600 .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
601 */
602 PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
603 {
604   PetscErrorCode           ierr;
605   PetscInt                 offsetp;
606   DMNetworkComponentHeader header;
607   DM_Network               *network = (DM_Network*)dm->data;
608 
609   PetscFunctionBegin;
610   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
611   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
612   if (compkey) *compkey = header->key[compnum];
613   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
614   PetscFunctionReturn(0);
615 }
616 
617 /*@
618   DMNetworkGetComponent - Returns the network component and its key
619 
620   Not Collective
621 
622   Input Parameters
623 + dm - DMNetwork object
624 . p  - edge or vertex point
625 - compnum - component number
626 
627   Output Parameters:
628 + compkey - the key set for this computing during registration
629 - component - the component data
630 
631   Notes:
632   Typical usage:
633 
634   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
635   Loop over vertices or edges
636     DMNetworkGetNumComponents(dm,v,&numcomps);
637     Loop over numcomps
638       DMNetworkGetComponent(dm,v,compnum,&key,&component);
639 
640   Level: intermediate
641 
642 .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset
643 @*/
644 PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component)
645 {
646   PetscErrorCode ierr;
647   DM_Network     *network = (DM_Network*)dm->data;
648   PetscInt       offsetd;
649 
650   PetscFunctionBegin;
651 
652   ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr);
653   *component = network->componentdataarray+offsetd;
654 
655   PetscFunctionReturn(0);
656 }
657 
658 /*@
659   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
660 
661   Not Collective
662 
663   Input Parameters:
664 + dm           - The DMNetwork object
665 . p            - vertex/edge point
666 . componentkey - component key returned while registering the component
667 - compvalue    - pointer to the data structure for the component
668 
669   Level: intermediate
670 
671 .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
672 @*/
673 PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
674 {
675   DM_Network               *network = (DM_Network*)dm->data;
676   DMNetworkComponent       *component = &network->component[componentkey];
677   DMNetworkComponentHeader header = &network->header[p];
678   DMNetworkComponentValue  cvalue = &network->cvalue[p];
679   PetscErrorCode           ierr;
680 
681   PetscFunctionBegin;
682   if (header->ndata == MAX_DATA_AT_POINT) SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components at a point exceeds the max %D",MAX_DATA_AT_POINT);
683 
684   header->size[header->ndata] = component->size;
685   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
686   header->key[header->ndata] = componentkey;
687   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
688 
689   cvalue->data[header->ndata] = (void*)compvalue;
690   header->ndata++;
691   PetscFunctionReturn(0);
692 }
693 
694 /*@
695   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
696 
697   Not Collective
698 
699   Input Parameters:
700 + dm - The DMNetwork object
701 . p  - vertex/edge point
702 
703   Output Parameters:
704 . numcomponents - Number of components at the vertex/edge
705 
706   Level: intermediate
707 
708 .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
709 @*/
710 PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
711 {
712   PetscErrorCode ierr;
713   PetscInt       offset;
714   DM_Network     *network = (DM_Network*)dm->data;
715 
716   PetscFunctionBegin;
717   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
718   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
719   PetscFunctionReturn(0);
720 }
721 
722 /*@
723   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
724 
725   Not Collective
726 
727   Input Parameters:
728 + dm     - The DMNetwork object
729 - p      - the edge/vertex point
730 
731   Output Parameters:
732 . offset - the offset
733 
734   Level: intermediate
735 
736 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
737 @*/
738 PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
739 {
740   PetscErrorCode ierr;
741   DM_Network     *network = (DM_Network*)dm->data;
742 
743   PetscFunctionBegin;
744   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
745   PetscFunctionReturn(0);
746 }
747 
748 /*@
749   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
750 
751   Not Collective
752 
753   Input Parameters:
754 + dm      - The DMNetwork object
755 - p       - the edge/vertex point
756 
757   Output Parameters:
758 . offsetg - the offset
759 
760   Level: intermediate
761 
762 .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
763 @*/
764 PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
765 {
766   PetscErrorCode ierr;
767   DM_Network     *network = (DM_Network*)dm->data;
768 
769   PetscFunctionBegin;
770   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
771   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
772   PetscFunctionReturn(0);
773 }
774 
775 /*@
776   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
777 
778   Not Collective
779 
780   Input Parameters:
781 + dm     - The DMNetwork object
782 - p      - the edge point
783 
784   Output Parameters:
785 . offset - the offset
786 
787   Level: intermediate
788 
789 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
790 @*/
791 PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
792 {
793   PetscErrorCode ierr;
794   DM_Network     *network = (DM_Network*)dm->data;
795 
796   PetscFunctionBegin;
797 
798   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
799   PetscFunctionReturn(0);
800 }
801 
802 /*@
803   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
804 
805   Not Collective
806 
807   Input Parameters:
808 + dm     - The DMNetwork object
809 - p      - the vertex point
810 
811   Output Parameters:
812 . offset - the offset
813 
814   Level: intermediate
815 
816 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
817 @*/
818 PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
819 {
820   PetscErrorCode ierr;
821   DM_Network     *network = (DM_Network*)dm->data;
822 
823   PetscFunctionBegin;
824 
825   p -= network->vStart;
826 
827   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
828   PetscFunctionReturn(0);
829 }
830 /*@
831   DMNetworkAddNumVariables - Add number of variables associated with a given point.
832 
833   Not Collective
834 
835   Input Parameters:
836 + dm   - The DMNetworkObject
837 . p    - the vertex/edge point
838 - nvar - number of additional variables
839 
840   Level: intermediate
841 
842 .seealso: DMNetworkSetNumVariables
843 @*/
844 PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
845 {
846   PetscErrorCode ierr;
847   DM_Network     *network = (DM_Network*)dm->data;
848 
849   PetscFunctionBegin;
850   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
851   PetscFunctionReturn(0);
852 }
853 
854 /*@
855   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
856 
857   Not Collective
858 
859   Input Parameters:
860 + dm   - The DMNetworkObject
861 - p    - the vertex/edge point
862 
863   Output Parameters:
864 . nvar - number of variables
865 
866   Level: intermediate
867 
868 .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
869 @*/
870 PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
871 {
872   PetscErrorCode ierr;
873   DM_Network     *network = (DM_Network*)dm->data;
874 
875   PetscFunctionBegin;
876   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
877   PetscFunctionReturn(0);
878 }
879 
880 /*@
881   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
882 
883   Not Collective
884 
885   Input Parameters:
886 + dm   - The DMNetworkObject
887 . p    - the vertex/edge point
888 - nvar - number of variables
889 
890   Level: intermediate
891 
892 .seealso: DMNetworkAddNumVariables
893 @*/
894 PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
895 {
896   PetscErrorCode ierr;
897   DM_Network     *network = (DM_Network*)dm->data;
898 
899   PetscFunctionBegin;
900   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
901   PetscFunctionReturn(0);
902 }
903 
904 /* Sets up the array that holds the data for all components and its associated section. This
905    function is called during DMSetUp() */
906 PetscErrorCode DMNetworkComponentSetUp(DM dm)
907 {
908   PetscErrorCode              ierr;
909   DM_Network     *network = (DM_Network*)dm->data;
910   PetscInt                    arr_size;
911   PetscInt                    p,offset,offsetp;
912   DMNetworkComponentHeader header;
913   DMNetworkComponentValue  cvalue;
914   DMNetworkComponentGenericDataType      *componentdataarray;
915   PetscInt ncomp, i;
916 
917   PetscFunctionBegin;
918   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
919   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
920   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
921   componentdataarray = network->componentdataarray;
922   for (p = network->pStart; p < network->pEnd; p++) {
923     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
924     /* Copy header */
925     header = &network->header[p];
926     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
927     /* Copy data */
928     cvalue = &network->cvalue[p];
929     ncomp = header->ndata;
930     for (i = 0; i < ncomp; i++) {
931       offset = offsetp + network->dataheadersize + header->offset[i];
932       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
933     }
934   }
935   PetscFunctionReturn(0);
936 }
937 
938 /* Sets up the section for dofs. This routine is called during DMSetUp() */
939 PetscErrorCode DMNetworkVariablesSetUp(DM dm)
940 {
941   PetscErrorCode ierr;
942   DM_Network     *network = (DM_Network*)dm->data;
943 
944   PetscFunctionBegin;
945   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
946   PetscFunctionReturn(0);
947 }
948 
949 /*@C
950   DMNetworkGetComponentDataArray - Returns the component data array
951 
952   Not Collective
953 
954   Input Parameters:
955 . dm - The DMNetwork Object
956 
957   Output Parameters:
958 . componentdataarray - array that holds data for all components
959 
960   Level: intermediate
961 
962 .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
963 @*/
964 PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
965 {
966   DM_Network     *network = (DM_Network*)dm->data;
967 
968   PetscFunctionBegin;
969   *componentdataarray = network->componentdataarray;
970   PetscFunctionReturn(0);
971 }
972 
973 /* Get a subsection from a range of points */
974 PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
975 {
976   PetscErrorCode ierr;
977   PetscInt       i, nvar;
978 
979   PetscFunctionBegin;
980   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
981   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
982   for (i = pstart; i < pend; i++) {
983     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
984     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
985   }
986 
987   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
988   PetscFunctionReturn(0);
989 }
990 
991 /* Create a submap of points with a GlobalToLocal structure */
992 PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
993 {
994   PetscErrorCode ierr;
995   PetscInt       i, *subpoints;
996 
997   PetscFunctionBegin;
998   /* Create index sets to map from "points" to "subpoints" */
999   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
1000   for (i = pstart; i < pend; i++) {
1001     subpoints[i - pstart] = i;
1002   }
1003   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
1004   ierr = PetscFree(subpoints);CHKERRQ(ierr);
1005   PetscFunctionReturn(0);
1006 }
1007 
1008 /*@
1009   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
1010 
1011   Collective
1012 
1013   Input Parameters:
1014 . dm   - The DMNetworkObject
1015 
1016   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
1017 
1018   points = [0 1 2 3 4 5 6]
1019 
1020   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
1021 
1022   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
1023 
1024   Level: intermediate
1025 
1026 @*/
1027 PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
1028 {
1029   PetscErrorCode ierr;
1030   MPI_Comm       comm;
1031   PetscMPIInt    rank, size;
1032   DM_Network     *network = (DM_Network*)dm->data;
1033 
1034   PetscFunctionBegin;
1035   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1036   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
1037   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1038 
1039   /* Create maps for vertices and edges */
1040   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
1041   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
1042 
1043   /* Create local sub-sections */
1044   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
1045   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
1046 
1047   if (size > 1) {
1048     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
1049     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
1050   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
1051   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
1052   } else {
1053   /* create structures for vertex */
1054   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
1055   /* create structures for edge */
1056   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
1057   }
1058 
1059 
1060   /* Add viewers */
1061   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
1062   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
1063   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
1064   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
1065 
1066   PetscFunctionReturn(0);
1067 }
1068 
1069 /*@
1070   DMNetworkDistribute - Distributes the network and moves associated component data.
1071 
1072   Collective
1073 
1074   Input Parameter:
1075 + DM - the DMNetwork object
1076 - overlap - The overlap of partitions, 0 is the default
1077 
1078   Notes:
1079   Distributes the network with <overlap>-overlapping partitioning of the edges.
1080 
1081   Level: intermediate
1082 
1083 .seealso: DMNetworkCreate
1084 @*/
1085 PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
1086 {
1087   MPI_Comm       comm;
1088   PetscErrorCode ierr;
1089   PetscMPIInt    size;
1090   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
1091   DM_Network     *newDMnetwork;
1092   PetscSF        pointsf;
1093   DM             newDM;
1094   PetscPartitioner part;
1095   PetscInt         j,e,v,offset;
1096   DMNetworkComponentHeader header;
1097 
1098   PetscFunctionBegin;
1099 
1100   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
1101   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
1102   if (size == 1) PetscFunctionReturn(0);
1103 
1104   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
1105   newDMnetwork = (DM_Network*)newDM->data;
1106   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
1107 
1108   /* Enable runtime options for petscpartitioner */
1109   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
1110   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
1111 
1112   /* Distribute plex dm and dof section */
1113   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
1114 
1115   /* Distribute dof section */
1116   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
1117   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
1118   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
1119 
1120   /* Distribute data and associated section */
1121   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
1122 
1123   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
1124   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
1125   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
1126   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
1127   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
1128   newDMnetwork->NVertices = oldDMnetwork->NVertices;
1129   newDMnetwork->NEdges = oldDMnetwork->NEdges;
1130 
1131   /* Set Dof section as the default section for dm */
1132   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
1133   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
1134 
1135   /* Set up subnetwork info in the newDM */
1136   newDMnetwork->nsubnet = oldDMnetwork->nsubnet;
1137   ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr);
1138   /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already
1139      calculated in DMNetworkLayoutSetUp()
1140   */
1141   for(j=0; j < newDMnetwork->nsubnet; j++) {
1142     newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx;
1143     newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge;
1144   }
1145 
1146   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1147     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1148     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1149     newDMnetwork->subnet[header->subnetid].nedge++;
1150   }
1151 
1152   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1153     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1154     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1155     newDMnetwork->subnet[header->subnetid].nvtx++;
1156   }
1157 
1158   /* Now create the vertices and edge arrays for the subnetworks */
1159   for(j=0; j < newDMnetwork->nsubnet; j++) {
1160     ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr);
1161     ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr);
1162     /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop.
1163        These get updated when the vertices and edges are added. */
1164     newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0;
1165   }
1166 
1167   /* Set the vertices and edges in each subnetwork */
1168   for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) {
1169     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr);
1170     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1171     newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++]  = e;
1172   }
1173 
1174   for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) {
1175     ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr);
1176     header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr);
1177     newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++]  = v;
1178   }
1179 
1180   /* Destroy point SF */
1181   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
1182 
1183   ierr = DMDestroy(dm);CHKERRQ(ierr);
1184   *dm  = newDM;
1185   PetscFunctionReturn(0);
1186 }
1187 
1188 /*@C
1189   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
1190 
1191   Input Parameters:
1192 + masterSF - the original SF structure
1193 - map      - a ISLocalToGlobal mapping that contains the subset of points
1194 
1195   Output Parameters:
1196 . subSF    - a subset of the masterSF for the desired subset.
1197 */
1198 
1199 PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
1200 
1201   PetscErrorCode        ierr;
1202   PetscInt              nroots, nleaves, *ilocal_sub;
1203   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
1204   PetscInt              *local_points, *remote_points;
1205   PetscSFNode           *iremote_sub;
1206   const PetscInt        *ilocal;
1207   const PetscSFNode     *iremote;
1208 
1209   PetscFunctionBegin;
1210   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
1211 
1212   /* Look for leaves that pertain to the subset of points. Get the local ordering */
1213   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
1214   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
1215   for (i = 0; i < nleaves; i++) {
1216     if (ilocal_map[i] != -1) nleaves_sub += 1;
1217   }
1218   /* Re-number ilocal with subset numbering. Need information from roots */
1219   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
1220   for (i = 0; i < nroots; i++) local_points[i] = i;
1221   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
1222   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
1223   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
1224   /* Fill up graph using local (that is, local to the subset) numbering. */
1225   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
1226   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
1227   nleaves_sub = 0;
1228   for (i = 0; i < nleaves; i++) {
1229     if (ilocal_map[i] != -1) {
1230       ilocal_sub[nleaves_sub] = ilocal_map[i];
1231       iremote_sub[nleaves_sub].rank = iremote[i].rank;
1232       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
1233       nleaves_sub += 1;
1234     }
1235   }
1236   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
1237   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
1238 
1239   /* Create new subSF */
1240   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
1241   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
1242   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
1243   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
1244   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
1245   PetscFunctionReturn(0);
1246 }
1247 
1248 /*@C
1249   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
1250 
1251   Not Collective
1252 
1253   Input Parameters:
1254 + dm - The DMNetwork object
1255 - p  - the vertex point
1256 
1257   Output Paramters:
1258 + nedges - number of edges connected to this vertex point
1259 - edges  - List of edge points
1260 
1261   Level: intermediate
1262 
1263   Fortran Notes:
1264   Since it returns an array, this routine is only available in Fortran 90, and you must
1265   include petsc.h90 in your code.
1266 
1267 .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
1268 @*/
1269 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
1270 {
1271   PetscErrorCode ierr;
1272   DM_Network     *network = (DM_Network*)dm->data;
1273 
1274   PetscFunctionBegin;
1275   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
1276   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
1277   PetscFunctionReturn(0);
1278 }
1279 
1280 /*@C
1281   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
1282 
1283   Not Collective
1284 
1285   Input Parameters:
1286 + dm - The DMNetwork object
1287 - p  - the edge point
1288 
1289   Output Paramters:
1290 . vertices  - vertices connected to this edge
1291 
1292   Level: intermediate
1293 
1294   Fortran Notes:
1295   Since it returns an array, this routine is only available in Fortran 90, and you must
1296   include petsc.h90 in your code.
1297 
1298 .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
1299 @*/
1300 PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
1301 {
1302   PetscErrorCode ierr;
1303   DM_Network     *network = (DM_Network*)dm->data;
1304 
1305   PetscFunctionBegin;
1306   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
1307   PetscFunctionReturn(0);
1308 }
1309 
1310 /*@
1311   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
1312 
1313   Not Collective
1314 
1315   Input Parameters:
1316 + dm - The DMNetwork object
1317 . p  - the vertex point
1318 
1319   Output Parameter:
1320 . isghost - TRUE if the vertex is a ghost point
1321 
1322   Level: intermediate
1323 
1324 .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
1325 @*/
1326 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
1327 {
1328   PetscErrorCode ierr;
1329   DM_Network     *network = (DM_Network*)dm->data;
1330   PetscInt       offsetg;
1331   PetscSection   sectiong;
1332 
1333   PetscFunctionBegin;
1334   *isghost = PETSC_FALSE;
1335   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
1336   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
1337   if (offsetg < 0) *isghost = PETSC_TRUE;
1338   PetscFunctionReturn(0);
1339 }
1340 
1341 PetscErrorCode DMSetUp_Network(DM dm)
1342 {
1343   PetscErrorCode ierr;
1344   DM_Network     *network=(DM_Network*)dm->data;
1345 
1346   PetscFunctionBegin;
1347   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
1348   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
1349 
1350   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
1351   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
1352   PetscFunctionReturn(0);
1353 }
1354 
1355 /*@
1356     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
1357                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
1358 
1359     Collective
1360 
1361     Input Parameters:
1362 +   dm - The DMNetwork object
1363 .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
1364 -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
1365 
1366     Level: intermediate
1367 
1368 @*/
1369 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
1370 {
1371   DM_Network     *network=(DM_Network*)dm->data;
1372   PetscErrorCode ierr;
1373 
1374   PetscFunctionBegin;
1375   network->userEdgeJacobian   = eflg;
1376   network->userVertexJacobian = vflg;
1377 
1378   if (eflg && !network->Je) {
1379     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
1380   }
1381 
1382   if (vflg && !network->Jv) {
1383     PetscInt       i,*vptr,nedges,vStart=network->vStart;
1384     PetscInt       nVertices = network->nVertices,nedges_total;
1385     const PetscInt *edges;
1386 
1387     /* count nvertex_total */
1388     nedges_total = 0;
1389     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
1390 
1391     vptr[0] = 0;
1392     for (i=0; i<nVertices; i++) {
1393       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1394       nedges_total += nedges;
1395       vptr[i+1] = vptr[i] + 2*nedges + 1;
1396     }
1397 
1398     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
1399     network->Jvptr = vptr;
1400   }
1401   PetscFunctionReturn(0);
1402 }
1403 
1404 /*@
1405     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
1406 
1407     Not Collective
1408 
1409     Input Parameters:
1410 +   dm - The DMNetwork object
1411 .   p  - the edge point
1412 -   J - array (size = 3) of Jacobian submatrices for this edge point:
1413         J[0]: this edge
1414         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
1415 
1416     Level: intermediate
1417 
1418 .seealso: DMNetworkVertexSetMatrix
1419 @*/
1420 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
1421 {
1422   DM_Network     *network=(DM_Network*)dm->data;
1423 
1424   PetscFunctionBegin;
1425   if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
1426 
1427   if (J) {
1428     network->Je[3*p]   = J[0];
1429     network->Je[3*p+1] = J[1];
1430     network->Je[3*p+2] = J[2];
1431   }
1432   PetscFunctionReturn(0);
1433 }
1434 
1435 /*@
1436     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
1437 
1438     Not Collective
1439 
1440     Input Parameters:
1441 +   dm - The DMNetwork object
1442 .   p  - the vertex point
1443 -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
1444         J[0]:       this vertex
1445         J[1+2*i]:   i-th supporting edge
1446         J[1+2*i+1]: i-th connected vertex
1447 
1448     Level: intermediate
1449 
1450 .seealso: DMNetworkEdgeSetMatrix
1451 @*/
1452 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
1453 {
1454   PetscErrorCode ierr;
1455   DM_Network     *network=(DM_Network*)dm->data;
1456   PetscInt       i,*vptr,nedges,vStart=network->vStart;
1457   const PetscInt *edges;
1458 
1459   PetscFunctionBegin;
1460   if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1461 
1462   if (J) {
1463     vptr = network->Jvptr;
1464     network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
1465 
1466     /* Set Jacobian for each supporting edge and connected vertex */
1467     ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1468     for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1469   }
1470   PetscFunctionReturn(0);
1471 }
1472 
1473 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1474 {
1475   PetscErrorCode ierr;
1476   PetscInt       j;
1477   PetscScalar    val=(PetscScalar)ncols;
1478 
1479   PetscFunctionBegin;
1480   if (!ghost) {
1481     for (j=0; j<nrows; j++) {
1482       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1483     }
1484   } else {
1485     for (j=0; j<nrows; j++) {
1486       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1487     }
1488   }
1489   PetscFunctionReturn(0);
1490 }
1491 
1492 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1493 {
1494   PetscErrorCode ierr;
1495   PetscInt       j,ncols_u;
1496   PetscScalar    val;
1497 
1498   PetscFunctionBegin;
1499   if (!ghost) {
1500     for (j=0; j<nrows; j++) {
1501       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1502       val = (PetscScalar)ncols_u;
1503       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1504       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1505     }
1506   } else {
1507     for (j=0; j<nrows; j++) {
1508       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1509       val = (PetscScalar)ncols_u;
1510       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1511       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1512     }
1513   }
1514   PetscFunctionReturn(0);
1515 }
1516 
1517 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1518 {
1519   PetscErrorCode ierr;
1520 
1521   PetscFunctionBegin;
1522   if (Ju) {
1523     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
1524   } else {
1525     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
1526   }
1527   PetscFunctionReturn(0);
1528 }
1529 
1530 PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1531 {
1532   PetscErrorCode ierr;
1533   PetscInt       j,*cols;
1534   PetscScalar    *zeros;
1535 
1536   PetscFunctionBegin;
1537   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1538   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1539   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1540   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
1541   PetscFunctionReturn(0);
1542 }
1543 
1544 PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1545 {
1546   PetscErrorCode ierr;
1547   PetscInt       j,M,N,row,col,ncols_u;
1548   const PetscInt *cols;
1549   PetscScalar    zero=0.0;
1550 
1551   PetscFunctionBegin;
1552   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
1553   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);
1554 
1555   for (row=0; row<nrows; row++) {
1556     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
1557     for (j=0; j<ncols_u; j++) {
1558       col = cols[j] + cstart;
1559       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
1560     }
1561     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
1562   }
1563   PetscFunctionReturn(0);
1564 }
1565 
1566 PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1567 {
1568   PetscErrorCode ierr;
1569 
1570   PetscFunctionBegin;
1571   if (Ju) {
1572     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1573   } else {
1574     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1575   }
1576   PetscFunctionReturn(0);
1577 }
1578 
1579 /* 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.
1580 */
1581 PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
1582 {
1583   PetscErrorCode ierr;
1584   PetscInt       i, size, dof;
1585   PetscInt       *glob2loc;
1586 
1587   PetscFunctionBegin;
1588   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
1589   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
1590 
1591   for (i = 0; i < size; i++) {
1592     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
1593     dof = (dof >= 0) ? dof : -(dof + 1);
1594     glob2loc[i] = dof;
1595   }
1596 
1597   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
1598 #if 0
1599   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1600 #endif
1601   PetscFunctionReturn(0);
1602 }
1603 
1604 #include <petsc/private/matimpl.h>
1605 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
1606 {
1607   PetscErrorCode ierr;
1608   PetscMPIInt    rank, size;
1609   DM_Network     *network = (DM_Network*) dm->data;
1610   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1611   PetscInt       cstart,ncols,j,e,v;
1612   PetscBool      ghost,ghost_vc,ghost2,isNest;
1613   Mat            Juser;
1614   PetscSection   sectionGlobal;
1615   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1616   const PetscInt *edges,*cone;
1617   MPI_Comm       comm;
1618   MatType        mtype;
1619   Vec            vd_nz,vo_nz;
1620   PetscInt       *dnnz,*onnz;
1621   PetscScalar    *vdnz,*vonz;
1622 
1623   PetscFunctionBegin;
1624   mtype = dm->mattype;
1625   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
1626 
1627   if (isNest) {
1628     /* ierr = DMCreateMatrix_Network_Nest(); */
1629     PetscInt   eDof, vDof;
1630     Mat        j11, j12, j21, j22, bA[2][2];
1631     ISLocalToGlobalMapping eISMap, vISMap;
1632 
1633     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1634     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1635     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1636 
1637     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
1638     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
1639 
1640     ierr = MatCreate(comm, &j11);CHKERRQ(ierr);
1641     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1642     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
1643 
1644     ierr = MatCreate(comm, &j12);CHKERRQ(ierr);
1645     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
1646     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
1647 
1648     ierr = MatCreate(comm, &j21);CHKERRQ(ierr);
1649     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1650     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
1651 
1652     ierr = MatCreate(comm, &j22);CHKERRQ(ierr);
1653     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1654     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
1655 
1656     bA[0][0] = j11;
1657     bA[0][1] = j12;
1658     bA[1][0] = j21;
1659     bA[1][1] = j22;
1660 
1661     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
1662     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
1663 
1664     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
1665     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
1666     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
1667     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
1668 
1669     ierr = MatSetUp(j11);CHKERRQ(ierr);
1670     ierr = MatSetUp(j12);CHKERRQ(ierr);
1671     ierr = MatSetUp(j21);CHKERRQ(ierr);
1672     ierr = MatSetUp(j22);CHKERRQ(ierr);
1673 
1674     ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
1675     ierr = MatSetUp(*J);CHKERRQ(ierr);
1676     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
1677     ierr = MatDestroy(&j11);CHKERRQ(ierr);
1678     ierr = MatDestroy(&j12);CHKERRQ(ierr);
1679     ierr = MatDestroy(&j21);CHKERRQ(ierr);
1680     ierr = MatDestroy(&j22);CHKERRQ(ierr);
1681 
1682     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1683     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1684 
1685     /* Free structures */
1686     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
1687     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
1688 
1689     PetscFunctionReturn(0);
1690   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1691     /* user does not provide Jacobian blocks */
1692     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1693     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1694     PetscFunctionReturn(0);
1695   }
1696 
1697   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
1698   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1699   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1700   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1701 
1702   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
1703   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
1704 
1705   /* (1) Set matrix preallocation */
1706   /*------------------------------*/
1707   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1708   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1709   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1710   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1711   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1712   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1713 
1714   /* Set preallocation for edges */
1715   /*-----------------------------*/
1716   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1717 
1718   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1719   for (e=eStart; e<eEnd; e++) {
1720     /* Get row indices */
1721     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1722     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1723     if (nrows) {
1724       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1725 
1726       /* Set preallocation for conntected vertices */
1727       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1728       for (v=0; v<2; v++) {
1729         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1730 
1731         if (network->Je) {
1732           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1733         } else Juser = NULL;
1734         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
1735         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1736       }
1737 
1738       /* Set preallocation for edge self */
1739       cstart = rstart;
1740       if (network->Je) {
1741         Juser = network->Je[3*e]; /* Jacobian(e,e) */
1742       } else Juser = NULL;
1743       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1744     }
1745   }
1746 
1747   /* Set preallocation for vertices */
1748   /*--------------------------------*/
1749   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1750   if (vEnd - vStart) vptr = network->Jvptr;
1751 
1752   for (v=vStart; v<vEnd; v++) {
1753     /* Get row indices */
1754     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1755     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1756     if (!nrows) continue;
1757 
1758     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1759     if (ghost) {
1760       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1761     } else {
1762       rows_v = rows;
1763     }
1764 
1765     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1766 
1767     /* Get supporting edges and connected vertices */
1768     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1769 
1770     for (e=0; e<nedges; e++) {
1771       /* Supporting edges */
1772       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1773       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1774 
1775       if (network->Jv) {
1776         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1777       } else Juser = NULL;
1778       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1779 
1780       /* Connected vertices */
1781       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1782       vc = (v == cone[0]) ? cone[1]:cone[0];
1783       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1784 
1785       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1786 
1787       if (network->Jv) {
1788         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1789       } else Juser = NULL;
1790       if (ghost_vc||ghost) {
1791         ghost2 = PETSC_TRUE;
1792       } else {
1793         ghost2 = PETSC_FALSE;
1794       }
1795       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1796     }
1797 
1798     /* Set preallocation for vertex self */
1799     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1800     if (!ghost) {
1801       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1802       if (network->Jv) {
1803         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1804       } else Juser = NULL;
1805       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1806     }
1807     if (ghost) {
1808       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1809     }
1810   }
1811 
1812   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1813   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
1814 
1815   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
1816 
1817   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1818   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1819 
1820   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1821   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1822   for (j=0; j<localSize; j++) {
1823     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1824     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1825   }
1826   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1827   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1828   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1829   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1830 
1831   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
1832   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
1833   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1834 
1835   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
1836 
1837   /* (2) Set matrix entries for edges */
1838   /*----------------------------------*/
1839   for (e=eStart; e<eEnd; e++) {
1840     /* Get row indices */
1841     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1842     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1843     if (nrows) {
1844       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1845 
1846       /* Set matrix entries for conntected vertices */
1847       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1848       for (v=0; v<2; v++) {
1849         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1850         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1851 
1852         if (network->Je) {
1853           Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1854         } else Juser = NULL;
1855         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1856       }
1857 
1858       /* Set matrix entries for edge self */
1859       cstart = rstart;
1860       if (network->Je) {
1861         Juser = network->Je[3*e]; /* Jacobian(e,e) */
1862       } else Juser = NULL;
1863       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1864     }
1865   }
1866 
1867   /* Set matrix entries for vertices */
1868   /*---------------------------------*/
1869   for (v=vStart; v<vEnd; v++) {
1870     /* Get row indices */
1871     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1872     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1873     if (!nrows) continue;
1874 
1875     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1876     if (ghost) {
1877       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1878     } else {
1879       rows_v = rows;
1880     }
1881     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1882 
1883     /* Get supporting edges and connected vertices */
1884     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1885 
1886     for (e=0; e<nedges; e++) {
1887       /* Supporting edges */
1888       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1889       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1890 
1891       if (network->Jv) {
1892         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1893       } else Juser = NULL;
1894       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1895 
1896       /* Connected vertices */
1897       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1898       vc = (v == cone[0]) ? cone[1]:cone[0];
1899 
1900       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
1901       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1902 
1903       if (network->Jv) {
1904         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1905       } else Juser = NULL;
1906       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1907     }
1908 
1909     /* Set matrix entries for vertex self */
1910     if (!ghost) {
1911       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1912       if (network->Jv) {
1913         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1914       } else Juser = NULL;
1915       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1916     }
1917     if (ghost) {
1918       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1919     }
1920   }
1921   ierr = PetscFree(rows);CHKERRQ(ierr);
1922 
1923   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1924   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1925 
1926   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1927   PetscFunctionReturn(0);
1928 }
1929 
1930 PetscErrorCode DMDestroy_Network(DM dm)
1931 {
1932   PetscErrorCode ierr;
1933   DM_Network     *network = (DM_Network*) dm->data;
1934   PetscInt       j;
1935 
1936   PetscFunctionBegin;
1937   if (--network->refct > 0) PetscFunctionReturn(0);
1938   if (network->Je) {
1939     ierr = PetscFree(network->Je);CHKERRQ(ierr);
1940   }
1941   if (network->Jv) {
1942     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
1943     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
1944   }
1945 
1946   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
1947   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
1948   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
1949   if (network->vertex.sf) {
1950     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
1951   }
1952   /* edge */
1953   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
1954   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
1955   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
1956   if (network->edge.sf) {
1957     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
1958   }
1959   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
1960   network->edges = NULL;
1961   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
1962   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
1963 
1964   for(j=0; j < network->nsubnet; j++) {
1965     ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr);
1966     ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr);
1967   }
1968   ierr = PetscFree(network->subnet);CHKERRQ(ierr);
1969   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
1970   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
1971   ierr = PetscFree(network->header);CHKERRQ(ierr);
1972   ierr = PetscFree(network);CHKERRQ(ierr);
1973   PetscFunctionReturn(0);
1974 }
1975 
1976 PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
1977 {
1978   PetscErrorCode ierr;
1979   DM_Network     *network = (DM_Network*) dm->data;
1980 
1981   PetscFunctionBegin;
1982   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
1983   PetscFunctionReturn(0);
1984 }
1985 
1986 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
1987 {
1988   PetscErrorCode ierr;
1989   DM_Network     *network = (DM_Network*) dm->data;
1990 
1991   PetscFunctionBegin;
1992   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
1993   PetscFunctionReturn(0);
1994 }
1995 
1996 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
1997 {
1998   PetscErrorCode ierr;
1999   DM_Network     *network = (DM_Network*) dm->data;
2000 
2001   PetscFunctionBegin;
2002   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
2003   PetscFunctionReturn(0);
2004 }
2005 
2006 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
2007 {
2008   PetscErrorCode ierr;
2009   DM_Network     *network = (DM_Network*) dm->data;
2010 
2011   PetscFunctionBegin;
2012   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
2013   PetscFunctionReturn(0);
2014 }
2015 
2016 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
2017 {
2018   PetscErrorCode ierr;
2019   DM_Network     *network = (DM_Network*) dm->data;
2020 
2021   PetscFunctionBegin;
2022   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
2023   PetscFunctionReturn(0);
2024 }
2025