xref: /petsc/src/dm/impls/network/network.c (revision 97a9ccd30fd95168d2009ddcb3942b90612f60fa)
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 local and global vertices and edges.
29 
30   Collective on DM
31 
32   Input Parameters:
33 + dm - the dm object
34 . nV - number of local vertices
35 . nE - number of local edges
36 . NV - number of global vertices (or PETSC_DETERMINE)
37 - NE - number of global edges (or PETSC_DETERMINE)
38 
39    Notes
40    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
41 
42    You cannot change the sizes once they have been set
43 
44    Level: intermediate
45 
46 .seealso: DMNetworkCreate()
47 @*/
48 PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE)
49 {
50   PetscErrorCode ierr;
51   DM_Network     *network = (DM_Network*) dm->data;
52   PetscInt       a[2],b[2];
53 
54   PetscFunctionBegin;
55   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
56   if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4);
57   if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5);
58   if (NV > 0 && nV > NV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local vertex size %D cannot be larger than global vertex size %D",nV,NV);
59   if (NE > 0 && nE > NE) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local edge size %D cannot be larger than global edge size %D",nE,NE);
60   if ((network->nNodes >= 0 || network->NNodes >= 0) && (network->nNodes != nV || network->NNodes != NV)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset vertex sizes to %D local %D global after previously setting them to %D local %D global",nV,NV,network->nNodes,network->NNodes);
61   if ((network->nEdges >= 0 || network->NEdges >= 0) && (network->nEdges != nE || network->NEdges != NE)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset edge sizes to %D local %D global after previously setting them to %D local %D global",nE,NE,network->nEdges,network->NEdges);
62   if (NE < 0 || NV < 0) {
63     a[0] = nV; a[1] = nE;
64     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
65     NV = b[0]; NE = b[1];
66   }
67   network->nNodes = nV;
68   network->NNodes = NV;
69   network->nEdges = nE;
70   network->NEdges = NE;
71   PetscFunctionReturn(0);
72 }
73 
74 /*@
75   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
76 
77   Logically collective on DM
78 
79   Input Parameters:
80 . edges - list of edges
81 
82   Notes:
83   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
84   not be destroyed before the call to DMNetworkLayoutSetUp
85 
86   Level: intermediate
87 
88 .seealso: DMNetworkCreate, DMNetworkSetSizes
89 @*/
90 PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[])
91 {
92   DM_Network *network = (DM_Network*) dm->data;
93 
94   PetscFunctionBegin;
95   network->edges = edgelist;
96   PetscFunctionReturn(0);
97 }
98 
99 /*@
100   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
101 
102   Collective on DM
103 
104   Input Parameters
105 . DM - the dmnetwork object
106 
107   Notes:
108   This routine should be called after the network sizes and edgelists have been provided. It creates
109   the bare layout of the network and sets up the network to begin insertion of components.
110 
111   All the components should be registered before calling this routine.
112 
113   Level: intermediate
114 
115 .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
116 @*/
117 PetscErrorCode DMNetworkLayoutSetUp(DM dm)
118 {
119   PetscErrorCode ierr;
120   DM_Network     *network = (DM_Network*) dm->data;
121   PetscInt       dim = 1; /* One dimensional network */
122   PetscInt       numCorners=2;
123   PetscInt       spacedim=2;
124   double         *vertexcoords=NULL;
125   PetscInt       i;
126   PetscInt       ndata;
127 
128   PetscFunctionBegin;
129   if (network->nNodes) {
130     ierr = PetscCalloc1(numCorners*network->nNodes,&vertexcoords);CHKERRQ(ierr);
131   }
132   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
133   if (network->nNodes) {
134     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
135   }
136   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
137   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
138   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
139 
140   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
141   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
142   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
143   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
144 
145 
146 
147   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
148   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
149   for (i = network->pStart; i < network->pEnd; i++) {
150     network->header[i].ndata = 0;
151     ndata = network->header[i].ndata;
152     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
153     network->header[i].offset[ndata] = 0;
154   }
155   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
156   PetscFunctionReturn(0);
157 }
158 
159 /*@C
160   DMNetworkRegisterComponent - Registers the network component
161 
162   Logically collective on DM
163 
164   Input Parameters
165 + dm   - the network object
166 . name - the component name
167 - size - the storage size in bytes for this component data
168 
169    Output Parameters
170 .   key - an integer key that defines the component
171 
172    Notes
173    This routine should be called by all processors before calling DMNetworkLayoutSetup().
174 
175    Level: intermediate
176 
177 .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
178 @*/
179 PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
180 {
181   PetscErrorCode        ierr;
182   DM_Network            *network = (DM_Network*) dm->data;
183   DMNetworkComponent    *component=&network->component[network->ncomponent];
184   PetscBool             flg=PETSC_FALSE;
185   PetscInt              i;
186 
187   PetscFunctionBegin;
188 
189   for (i=0; i < network->ncomponent; i++) {
190     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
191     if (flg) {
192       *key = i;
193       PetscFunctionReturn(0);
194     }
195   }
196 
197   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
198   component->size = size/sizeof(DMNetworkComponentGenericDataType);
199   *key = network->ncomponent;
200   network->ncomponent++;
201   PetscFunctionReturn(0);
202 }
203 
204 /*@
205   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
206 
207   Not Collective
208 
209   Input Parameters:
210 + dm - The DMNetwork object
211 
212   Output Paramters:
213 + vStart - The first vertex point
214 - vEnd   - One beyond the last vertex point
215 
216   Level: intermediate
217 
218 .seealso: DMNetworkGetEdgeRange
219 @*/
220 PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
221 {
222   DM_Network     *network = (DM_Network*)dm->data;
223 
224   PetscFunctionBegin;
225   if (vStart) *vStart = network->vStart;
226   if (vEnd) *vEnd = network->vEnd;
227   PetscFunctionReturn(0);
228 }
229 
230 /*@
231   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
232 
233   Not Collective
234 
235   Input Parameters:
236 + dm - The DMNetwork object
237 
238   Output Paramters:
239 + eStart - The first edge point
240 - eEnd   - One beyond the last edge point
241 
242   Level: intermediate
243 
244 .seealso: DMNetworkGetVertexRange
245 @*/
246 PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
247 {
248   DM_Network     *network = (DM_Network*)dm->data;
249 
250   PetscFunctionBegin;
251   if (eStart) *eStart = network->eStart;
252   if (eEnd) *eEnd = network->eEnd;
253   PetscFunctionReturn(0);
254 }
255 
256 /*@
257   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
258 
259   Not Collective
260 
261   Input Parameters:
262 + dm           - The DMNetwork object
263 . p            - vertex/edge point
264 . componentkey - component key returned while registering the component
265 - compvalue    - pointer to the data structure for the component
266 
267   Level: intermediate
268 
269 .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
270 @*/
271 PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
272 {
273   DM_Network               *network = (DM_Network*)dm->data;
274   DMNetworkComponent       *component = &network->component[componentkey];
275   DMNetworkComponentHeader header = &network->header[p];
276   DMNetworkComponentValue  cvalue = &network->cvalue[p];
277   PetscErrorCode           ierr;
278 
279   PetscFunctionBegin;
280   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);
281 
282   header->size[header->ndata] = component->size;
283   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
284   header->key[header->ndata] = componentkey;
285   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
286 
287   cvalue->data[header->ndata] = (void*)compvalue;
288   header->ndata++;
289   PetscFunctionReturn(0);
290 }
291 
292 /*@
293   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
294 
295   Not Collective
296 
297   Input Parameters:
298 + dm - The DMNetwork object
299 . p  - vertex/edge point
300 
301   Output Parameters:
302 . numcomponents - Number of components at the vertex/edge
303 
304   Level: intermediate
305 
306 .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
307 @*/
308 PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
309 {
310   PetscErrorCode ierr;
311   PetscInt       offset;
312   DM_Network     *network = (DM_Network*)dm->data;
313 
314   PetscFunctionBegin;
315   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
316   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
317   PetscFunctionReturn(0);
318 }
319 
320 /*@
321   DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the
322                                     component value from the component data array
323 
324   Not Collective
325 
326   Input Parameters:
327 + dm      - The DMNetwork object
328 . p       - vertex/edge point
329 - compnum - component number
330 
331   Output Parameters:
332 + compkey - the key obtained when registering the component
333 - offset  - offset into the component data array associated with the vertex/edge point
334 
335   Notes:
336   Typical usage:
337 
338   DMNetworkGetComponentDataArray(dm, &arr);
339   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
340   Loop over vertices or edges
341     DMNetworkGetNumComponents(dm,v,&numcomps);
342     Loop over numcomps
343       DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset);
344       compdata = (UserCompDataType)(arr+offset);
345 
346   Level: intermediate
347 
348 .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
349 @*/
350 PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
351 {
352   PetscErrorCode           ierr;
353   PetscInt                 offsetp;
354   DMNetworkComponentHeader header;
355   DM_Network               *network = (DM_Network*)dm->data;
356 
357   PetscFunctionBegin;
358   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
359   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
360   if (compkey) *compkey = header->key[compnum];
361   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
362   PetscFunctionReturn(0);
363 }
364 
365 /*@
366   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
367 
368   Not Collective
369 
370   Input Parameters:
371 + dm     - The DMNetwork object
372 - p      - the edge/vertex point
373 
374   Output Parameters:
375 . offset - the offset
376 
377   Level: intermediate
378 
379 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
380 @*/
381 PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
382 {
383   PetscErrorCode ierr;
384   DM_Network     *network = (DM_Network*)dm->data;
385 
386   PetscFunctionBegin;
387   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
388   PetscFunctionReturn(0);
389 }
390 
391 /*@
392   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
393 
394   Not Collective
395 
396   Input Parameters:
397 + dm      - The DMNetwork object
398 - p       - the edge/vertex point
399 
400   Output Parameters:
401 . offsetg - the offset
402 
403   Level: intermediate
404 
405 .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
406 @*/
407 PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
408 {
409   PetscErrorCode ierr;
410   DM_Network     *network = (DM_Network*)dm->data;
411 
412   PetscFunctionBegin;
413   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
414   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost node */
415   PetscFunctionReturn(0);
416 }
417 
418 /*@
419   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
420 
421   Not Collective
422 
423   Input Parameters:
424 + dm     - The DMNetwork object
425 - p      - the edge point
426 
427   Output Parameters:
428 . offset - the offset
429 
430   Level: intermediate
431 
432 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
433 @*/
434 PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
435 {
436   PetscErrorCode ierr;
437   DM_Network     *network = (DM_Network*)dm->data;
438 
439   PetscFunctionBegin;
440 
441   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
442   PetscFunctionReturn(0);
443 }
444 
445 /*@
446   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
447 
448   Not Collective
449 
450   Input Parameters:
451 + dm     - The DMNetwork object
452 - p      - the vertex point
453 
454   Output Parameters:
455 . offset - the offset
456 
457   Level: intermediate
458 
459 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
460 @*/
461 PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
462 {
463   PetscErrorCode ierr;
464   DM_Network     *network = (DM_Network*)dm->data;
465 
466   PetscFunctionBegin;
467 
468   p -= network->vStart;
469 
470   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
471   PetscFunctionReturn(0);
472 }
473 /*@
474   DMNetworkAddNumVariables - Add number of variables associated with a given point.
475 
476   Not Collective
477 
478   Input Parameters:
479 + dm   - The DMNetworkObject
480 . p    - the vertex/edge point
481 - nvar - number of additional variables
482 
483   Level: intermediate
484 
485 .seealso: DMNetworkSetNumVariables
486 @*/
487 PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
488 {
489   PetscErrorCode ierr;
490   DM_Network     *network = (DM_Network*)dm->data;
491 
492   PetscFunctionBegin;
493   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
494   PetscFunctionReturn(0);
495 }
496 
497 /*@
498   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
499 
500   Not Collective
501 
502   Input Parameters:
503 + dm   - The DMNetworkObject
504 - p    - the vertex/edge point
505 
506   Output Parameters:
507 . nvar - number of variables
508 
509   Level: intermediate
510 
511 .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
512 @*/
513 PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
514 {
515   PetscErrorCode ierr;
516   DM_Network     *network = (DM_Network*)dm->data;
517 
518   PetscFunctionBegin;
519   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
520   PetscFunctionReturn(0);
521 }
522 
523 /*@
524   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
525 
526   Not Collective
527 
528   Input Parameters:
529 + dm   - The DMNetworkObject
530 . p    - the vertex/edge point
531 - nvar - number of variables
532 
533   Level: intermediate
534 
535 .seealso: DMNetworkAddNumVariables
536 @*/
537 PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
538 {
539   PetscErrorCode ierr;
540   DM_Network     *network = (DM_Network*)dm->data;
541 
542   PetscFunctionBegin;
543   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
544   PetscFunctionReturn(0);
545 }
546 
547 /* Sets up the array that holds the data for all components and its associated section. This
548    function is called during DMSetUp() */
549 PetscErrorCode DMNetworkComponentSetUp(DM dm)
550 {
551   PetscErrorCode              ierr;
552   DM_Network     *network = (DM_Network*)dm->data;
553   PetscInt                    arr_size;
554   PetscInt                    p,offset,offsetp;
555   DMNetworkComponentHeader header;
556   DMNetworkComponentValue  cvalue;
557   DMNetworkComponentGenericDataType      *componentdataarray;
558   PetscInt ncomp, i;
559 
560   PetscFunctionBegin;
561   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
562   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
563   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
564   componentdataarray = network->componentdataarray;
565   for (p = network->pStart; p < network->pEnd; p++) {
566     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
567     /* Copy header */
568     header = &network->header[p];
569     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
570     /* Copy data */
571     cvalue = &network->cvalue[p];
572     ncomp = header->ndata;
573     for (i = 0; i < ncomp; i++) {
574       offset = offsetp + network->dataheadersize + header->offset[i];
575       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
576     }
577   }
578   PetscFunctionReturn(0);
579 }
580 
581 /* Sets up the section for dofs. This routine is called during DMSetUp() */
582 PetscErrorCode DMNetworkVariablesSetUp(DM dm)
583 {
584   PetscErrorCode ierr;
585   DM_Network     *network = (DM_Network*)dm->data;
586 
587   PetscFunctionBegin;
588   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
589   PetscFunctionReturn(0);
590 }
591 
592 /*@C
593   DMNetworkGetComponentDataArray - Returns the component data array
594 
595   Not Collective
596 
597   Input Parameters:
598 . dm - The DMNetwork Object
599 
600   Output Parameters:
601 . componentdataarray - array that holds data for all components
602 
603   Level: intermediate
604 
605 .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents
606 @*/
607 PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
608 {
609   DM_Network     *network = (DM_Network*)dm->data;
610 
611   PetscFunctionBegin;
612   *componentdataarray = network->componentdataarray;
613   PetscFunctionReturn(0);
614 }
615 
616 /* Get a subsection from a range of points */
617 PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
618 {
619   PetscErrorCode ierr;
620   PetscInt       i, nvar;
621 
622   PetscFunctionBegin;
623   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
624   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
625   for (i = pstart; i < pend; i++) {
626     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
627     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
628   }
629 
630   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
631   PetscFunctionReturn(0);
632 }
633 
634 /* Create a submap of points with a GlobalToLocal structure */
635 PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
636 {
637   PetscErrorCode ierr;
638   PetscInt       i, *subpoints;
639 
640   PetscFunctionBegin;
641   /* Create index sets to map from "points" to "subpoints" */
642   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
643   for (i = pstart; i < pend; i++) {
644     subpoints[i - pstart] = i;
645   }
646   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
647   ierr = PetscFree(subpoints);CHKERRQ(ierr);
648   PetscFunctionReturn(0);
649 }
650 
651 /*@
652   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
653 
654   Collective
655 
656   Input Parameters:
657 . dm   - The DMNetworkObject
658 
659   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
660 
661   points = [0 1 2 3 4 5 6]
662 
663   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
664 
665   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
666 
667   Level: intermediate
668 
669 @*/
670 PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
671 {
672   PetscErrorCode ierr;
673   MPI_Comm       comm;
674   PetscMPIInt    rank, size;
675   DM_Network     *network = (DM_Network*)dm->data;
676 
677   PetscFunctionBegin;
678   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
679   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
680   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
681 
682   /* Create maps for vertices and edges */
683   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
684   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
685 
686   /* Create local sub-sections */
687   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
688   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
689 
690   if (size > 1) {
691     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
692     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
693   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
694   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
695   } else {
696   /* create structures for vertex */
697   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
698   /* create structures for edge */
699   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
700   }
701 
702 
703   /* Add viewers */
704   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
705   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
706   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
707   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
708 
709   PetscFunctionReturn(0);
710 }
711 /*@
712   DMNetworkDistribute - Distributes the network and moves associated component data.
713 
714   Collective
715 
716   Input Parameter:
717 + oldDM - the original DMNetwork object
718 - overlap - The overlap of partitions, 0 is the default
719 
720   Output Parameter:
721 . distDM - the distributed DMNetwork object
722 
723   Notes:
724   This routine should be called only when using multiple processors.
725 
726   Distributes the network with <overlap>-overlapping partitioning of the edges.
727 
728   Level: intermediate
729 
730 .seealso: DMNetworkCreate
731 @*/
732 PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM)
733 {
734   PetscErrorCode ierr;
735   DM_Network     *oldDMnetwork = (DM_Network*)oldDM->data;
736   PetscSF        pointsf;
737   DM             newDM;
738   DM_Network     *newDMnetwork;
739 
740   PetscFunctionBegin;
741   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr);
742   newDMnetwork = (DM_Network*)newDM->data;
743   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
744   /* Distribute plex dm and dof section */
745   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
746   /* Distribute dof section */
747   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr);
748   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
749   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr);
750   /* Distribute data and associated section */
751   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
752 
753 
754   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
755   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
756   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
757   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
758   newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart;
759   newDMnetwork->NNodes = oldDMnetwork->NNodes;
760   newDMnetwork->NEdges = oldDMnetwork->NEdges;
761 
762   /* Set Dof section as the default section for dm */
763   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
764   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
765 
766   /* Destroy point SF */
767   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
768 
769   *distDM = newDM;
770   PetscFunctionReturn(0);
771 }
772 
773 /*@C
774   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
775 
776   Input Parameters:
777 + masterSF - the original SF structure
778 - map      - a ISLocalToGlobal mapping that contains the subset of points
779 
780   Output Parameters:
781 . subSF    - a subset of the masterSF for the desired subset.
782 */
783 
784 PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
785 
786   PetscErrorCode        ierr;
787   PetscInt              nroots, nleaves, *ilocal_sub;
788   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
789   PetscInt              *local_points, *remote_points;
790   PetscSFNode           *iremote_sub;
791   const PetscInt        *ilocal;
792   const PetscSFNode     *iremote;
793 
794   PetscFunctionBegin;
795   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
796 
797   /* Look for leaves that pertain to the subset of points. Get the local ordering */
798   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
799   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
800   for (i = 0; i < nleaves; i++) {
801     if (ilocal_map[i] != -1) nleaves_sub += 1;
802   }
803   /* Re-number ilocal with subset numbering. Need information from roots */
804   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
805   for (i = 0; i < nroots; i++) local_points[i] = i;
806   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
807   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
808   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
809   /* Fill up graph using local (that is, local to the subset) numbering. */
810   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
811   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
812   nleaves_sub = 0;
813   for (i = 0; i < nleaves; i++) {
814     if (ilocal_map[i] != -1) {
815       ilocal_sub[nleaves_sub] = ilocal_map[i];
816       iremote_sub[nleaves_sub].rank = iremote[i].rank;
817       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
818       nleaves_sub += 1;
819     }
820   }
821   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
822   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
823 
824   /* Create new subSF */
825   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
826   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
827   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
828   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
829   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
830 
831   PetscFunctionReturn(0);
832 }
833 
834 
835 /*@C
836   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
837 
838   Not Collective
839 
840   Input Parameters:
841 + dm - The DMNetwork object
842 - p  - the vertex point
843 
844   Output Paramters:
845 + nedges - number of edges connected to this vertex point
846 - edges  - List of edge points
847 
848   Level: intermediate
849 
850   Fortran Notes:
851   Since it returns an array, this routine is only available in Fortran 90, and you must
852   include petsc.h90 in your code.
853 
854 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
855 @*/
856 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
857 {
858   PetscErrorCode ierr;
859   DM_Network     *network = (DM_Network*)dm->data;
860 
861   PetscFunctionBegin;
862   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
863   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
864   PetscFunctionReturn(0);
865 }
866 
867 /*@C
868   DMNetworkGetConnectedNodes - Return the connected vertices for this edge point
869 
870   Not Collective
871 
872   Input Parameters:
873 + dm - The DMNetwork object
874 - p  - the edge point
875 
876   Output Paramters:
877 . vertices  - vertices connected to this edge
878 
879   Level: intermediate
880 
881   Fortran Notes:
882   Since it returns an array, this routine is only available in Fortran 90, and you must
883   include petsc.h90 in your code.
884 
885 .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
886 @*/
887 PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
888 {
889   PetscErrorCode ierr;
890   DM_Network     *network = (DM_Network*)dm->data;
891 
892   PetscFunctionBegin;
893   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
894   PetscFunctionReturn(0);
895 }
896 
897 /*@
898   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
899 
900   Not Collective
901 
902   Input Parameters:
903 + dm - The DMNetwork object
904 . p  - the vertex point
905 
906   Output Parameter:
907 . isghost - TRUE if the vertex is a ghost point
908 
909   Level: intermediate
910 
911 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
912 @*/
913 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
914 {
915   PetscErrorCode ierr;
916   DM_Network     *network = (DM_Network*)dm->data;
917   PetscInt       offsetg;
918   PetscSection   sectiong;
919 
920   PetscFunctionBegin;
921   *isghost = PETSC_FALSE;
922   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
923   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
924   if (offsetg < 0) *isghost = PETSC_TRUE;
925   PetscFunctionReturn(0);
926 }
927 
928 PetscErrorCode DMSetUp_Network(DM dm)
929 {
930   PetscErrorCode ierr;
931   DM_Network     *network=(DM_Network*)dm->data;
932 
933   PetscFunctionBegin;
934   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
935   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
936 
937   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
938   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
939   PetscFunctionReturn(0);
940 }
941 
942 /*@
943     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
944                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
945 
946     Collective
947 
948     Input Parameters:
949 +   dm - The DMNetwork object
950 .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
951 -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
952 
953     Level: intermediate
954 
955 @*/
956 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
957 {
958   DM_Network     *network=(DM_Network*)dm->data;
959 
960   PetscFunctionBegin;
961   network->userEdgeJacobian   = eflg;
962   network->userVertexJacobian = vflg;
963   PetscFunctionReturn(0);
964 }
965 
966 /*@
967     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
968 
969     Not Collective
970 
971     Input Parameters:
972 +   dm - The DMNetwork object
973 .   p  - the edge point
974 -   J - array (size = 3) of Jacobian submatrices for this edge point:
975         J[0]: this edge
976         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes()
977 
978     Level: intermediate
979 
980 .seealso: DMNetworkVertexSetMatrix
981 @*/
982 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
983 {
984   PetscErrorCode ierr;
985   DM_Network     *network=(DM_Network*)dm->data;
986 
987   PetscFunctionBegin;
988   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
989   if (!network->Je) {
990     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
991   }
992   network->Je[3*p]   = J[0];
993   network->Je[3*p+1] = J[1];
994   network->Je[3*p+2] = J[2];
995   PetscFunctionReturn(0);
996 }
997 
998 /*@
999     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
1000 
1001     Not Collective
1002 
1003     Input Parameters:
1004 +   dm - The DMNetwork object
1005 .   p  - the vertex point
1006 -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
1007         J[0]:       this vertex
1008         J[1+2*i]:   i-th supporting edge
1009         J[1+2*i+1]: i-th connected vertex
1010 
1011     Level: intermediate
1012 
1013 .seealso: DMNetworkEdgeSetMatrix
1014 @*/
1015 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
1016 {
1017   PetscErrorCode ierr;
1018   DM_Network     *network=(DM_Network*)dm->data;
1019   PetscInt       i,*vptr,nedges,vStart,vEnd;
1020   const PetscInt *edges;
1021 
1022   PetscFunctionBegin;
1023   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1024 
1025   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1026 
1027   if (!network->Jv) {
1028     PetscInt nNodes = network->nNodes,nedges_total;
1029     /* count nvertex_total */
1030     nedges_total = 0;
1031     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1032     ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr);
1033 
1034     vptr[0] = 0;
1035     for (i=0; i<nNodes; i++) {
1036       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1037       nedges_total += nedges;
1038       vptr[i+1] = vptr[i] + 2*nedges + 1;
1039     }
1040 
1041     ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr);
1042     network->Jvptr = vptr;
1043   }
1044 
1045   vptr = network->Jvptr;
1046   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
1047 
1048   /* Set Jacobian for each supporting edge and connected vertex */
1049   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1050   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1051   PetscFunctionReturn(0);
1052 }
1053 
1054 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1055 {
1056   PetscErrorCode ierr;
1057   PetscInt       j;
1058   PetscScalar    val=(PetscScalar)ncols;
1059 
1060   PetscFunctionBegin;
1061   if (!ghost) {
1062     for (j=0; j<nrows; j++) {
1063       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1064     }
1065   } else {
1066     for (j=0; j<nrows; j++) {
1067       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1068     }
1069   }
1070   PetscFunctionReturn(0);
1071 }
1072 
1073 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1074 {
1075   PetscErrorCode ierr;
1076   PetscInt       j,ncols_u;
1077   PetscScalar    val;
1078 
1079   PetscFunctionBegin;
1080   if (!ghost) {
1081     for (j=0; j<nrows; j++) {
1082       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1083       val = (PetscScalar)ncols_u;
1084       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1085       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1086     }
1087   } else {
1088     for (j=0; j<nrows; j++) {
1089       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1090       val = (PetscScalar)ncols_u;
1091       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1092       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1093     }
1094   }
1095   PetscFunctionReturn(0);
1096 }
1097 
1098 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1099 {
1100   PetscErrorCode ierr;
1101 
1102   PetscFunctionBegin;
1103   if (Ju) {
1104     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
1105   } else {
1106     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
1107   }
1108   PetscFunctionReturn(0);
1109 }
1110 
1111 PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1112 {
1113   PetscErrorCode ierr;
1114   PetscInt       j,*cols;
1115   PetscScalar    *zeros;
1116 
1117   PetscFunctionBegin;
1118   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1119   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1120   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1121   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
1122   PetscFunctionReturn(0);
1123 }
1124 
1125 PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1126 {
1127   PetscErrorCode ierr;
1128   PetscInt       j,M,N,row,col,ncols_u;
1129   const PetscInt *cols;
1130   PetscScalar    zero=0.0;
1131 
1132   PetscFunctionBegin;
1133   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
1134   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);
1135 
1136   for (row=0; row<nrows; row++) {
1137     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
1138     for (j=0; j<ncols_u; j++) {
1139       col = cols[j] + cstart;
1140       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
1141     }
1142     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
1143   }
1144   PetscFunctionReturn(0);
1145 }
1146 
1147 PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1148 {
1149   PetscErrorCode ierr;
1150 
1151   PetscFunctionBegin;
1152   if (Ju) {
1153     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1154   } else {
1155     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1156   }
1157   PetscFunctionReturn(0);
1158 }
1159 
1160 
1161 /* 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.
1162 */
1163 PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
1164 {
1165   PetscErrorCode ierr;
1166   PetscInt       i, size, dof;
1167   PetscInt       *glob2loc;
1168 
1169   PetscFunctionBegin;
1170   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
1171   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
1172 
1173   for (i = 0; i < size; i++) {
1174     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
1175     dof = (dof >= 0) ? dof : -(dof + 1);
1176     glob2loc[i] = dof;
1177   }
1178 
1179   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
1180 #if 0
1181   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1182 #endif
1183   PetscFunctionReturn(0);
1184 }
1185 
1186 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
1187 {
1188   PetscErrorCode ierr;
1189   PetscMPIInt    rank, size;
1190   DM_Network     *network = (DM_Network*) dm->data;
1191   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1192   PetscInt       cstart,ncols,j,e,v;
1193   PetscBool      ghost,ghost_vc,ghost2,isNest;
1194   Mat            Juser;
1195   PetscSection   sectionGlobal;
1196   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1197   const PetscInt *edges,*cone;
1198   MPI_Comm       comm;
1199   MatType        mtype;
1200   Vec            vd_nz,vo_nz;
1201   PetscInt       *dnnz,*onnz;
1202   PetscScalar    *vdnz,*vonz;
1203 
1204   PetscFunctionBegin;
1205   mtype = dm->mattype;
1206   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
1207 
1208   if (isNest) {
1209     /* ierr = DMCreateMatrix_Network_Nest(); */
1210     PetscInt   eDof, vDof;
1211     Mat        j11, j12, j21, j22, bA[2][2];
1212     ISLocalToGlobalMapping eISMap, vISMap;
1213 
1214     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1215     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1216     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1217 
1218     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
1219     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
1220 
1221     ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr);
1222     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1223     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
1224 
1225     ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr);
1226     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
1227     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
1228 
1229     ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr);
1230     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1231     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
1232 
1233     ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr);
1234     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1235     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
1236 
1237     bA[0][0] = j11;
1238     bA[0][1] = j12;
1239     bA[1][0] = j21;
1240     bA[1][1] = j22;
1241 
1242     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
1243     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
1244 
1245     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
1246     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
1247     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
1248     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
1249 
1250     ierr = MatSetUp(j11);CHKERRQ(ierr);
1251     ierr = MatSetUp(j12);CHKERRQ(ierr);
1252     ierr = MatSetUp(j21);CHKERRQ(ierr);
1253     ierr = MatSetUp(j22);CHKERRQ(ierr);
1254 
1255     ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
1256     ierr = MatSetUp(*J);CHKERRQ(ierr);
1257     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
1258     ierr = MatDestroy(&j11);CHKERRQ(ierr);
1259     ierr = MatDestroy(&j12);CHKERRQ(ierr);
1260     ierr = MatDestroy(&j21);CHKERRQ(ierr);
1261     ierr = MatDestroy(&j22);CHKERRQ(ierr);
1262 
1263     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1264     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1265 
1266     /* Free structures */
1267     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
1268     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
1269 
1270     PetscFunctionReturn(0);
1271   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1272     /* user does not provide Jacobian blocks */
1273     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1274     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1275     PetscFunctionReturn(0);
1276   }
1277 
1278   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
1279   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1280   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1281   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1282 
1283   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
1284   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
1285 
1286   /* (1) Set matrix preallocation */
1287   /*------------------------------*/
1288   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1289   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1290   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1291   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1292   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1293   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1294 
1295   /* Set preallocation for edges */
1296   /*-----------------------------*/
1297   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1298 
1299   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1300   for (e=eStart; e<eEnd; e++) {
1301     /* Get row indices */
1302     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1303     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1304     if (nrows) {
1305       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1306       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1307 
1308       /* Set preallocation for conntected vertices */
1309       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1310       for (v=0; v<2; v++) {
1311         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1312 
1313         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1314         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
1315         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1316       }
1317 
1318       /* Set preallocation for edge self */
1319       cstart = rstart;
1320       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1321       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1322     }
1323   }
1324 
1325   /* Set preallocation for vertices */
1326   /*--------------------------------*/
1327   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1328   if (vEnd - vStart) {
1329     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1330     vptr = network->Jvptr;
1331     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1332   }
1333 
1334   for (v=vStart; v<vEnd; v++) {
1335     /* Get row indices */
1336     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1337     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1338     if (!nrows) continue;
1339 
1340     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1341     if (ghost) {
1342       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1343     } else {
1344       rows_v = rows;
1345     }
1346 
1347     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1348 
1349     /* Get supporting edges and connected vertices */
1350     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1351 
1352     for (e=0; e<nedges; e++) {
1353       /* Supporting edges */
1354       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1355       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1356 
1357       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1358       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1359 
1360       /* Connected vertices */
1361       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1362       vc = (v == cone[0]) ? cone[1]:cone[0];
1363       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1364 
1365       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1366 
1367       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1368       if (ghost_vc||ghost) {
1369         ghost2 = PETSC_TRUE;
1370       } else {
1371         ghost2 = PETSC_FALSE;
1372       }
1373       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1374     }
1375 
1376     /* Set preallocation for vertex self */
1377     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1378     if (!ghost) {
1379       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1380       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1381       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1382     }
1383     if (ghost) {
1384       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1385     }
1386   }
1387 
1388   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1389   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
1390 
1391   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
1392 
1393   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1394   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1395 
1396   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1397   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1398   for (j=0; j<localSize; j++) {
1399     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1400     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1401   }
1402   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1403   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1404   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1405   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1406 
1407   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
1408   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
1409   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1410 
1411   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
1412 
1413   /* (2) Set matrix entries for edges */
1414   /*----------------------------------*/
1415   for (e=eStart; e<eEnd; e++) {
1416     /* Get row indices */
1417     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1418     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1419     if (nrows) {
1420       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1421 
1422       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1423 
1424       /* Set matrix entries for conntected vertices */
1425       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1426       for (v=0; v<2; v++) {
1427         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1428         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1429 
1430         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1431         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1432       }
1433 
1434       /* Set matrix entries for edge self */
1435       cstart = rstart;
1436       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1437       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1438     }
1439   }
1440 
1441   /* Set matrix entries for vertices */
1442   /*---------------------------------*/
1443   for (v=vStart; v<vEnd; v++) {
1444     /* Get row indices */
1445     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1446     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1447     if (!nrows) continue;
1448 
1449     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1450     if (ghost) {
1451       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1452     } else {
1453       rows_v = rows;
1454     }
1455     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1456 
1457     /* Get supporting edges and connected vertices */
1458     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1459 
1460     for (e=0; e<nedges; e++) {
1461       /* Supporting edges */
1462       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1463       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1464 
1465       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1466       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1467 
1468       /* Connected vertices */
1469       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1470       vc = (v == cone[0]) ? cone[1]:cone[0];
1471 
1472       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
1473       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1474 
1475       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1476       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1477     }
1478 
1479     /* Set matrix entries for vertex self */
1480     if (!ghost) {
1481       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1482       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1483       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1484     }
1485     if (ghost) {
1486       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1487     }
1488   }
1489   ierr = PetscFree(rows);CHKERRQ(ierr);
1490 
1491   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1492   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1493 
1494   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1495   PetscFunctionReturn(0);
1496 }
1497 
1498 PetscErrorCode DMDestroy_Network(DM dm)
1499 {
1500   PetscErrorCode ierr;
1501   DM_Network     *network = (DM_Network*) dm->data;
1502 
1503   PetscFunctionBegin;
1504   if (--network->refct > 0) PetscFunctionReturn(0);
1505   if (network->Je) {
1506     ierr = PetscFree(network->Je);CHKERRQ(ierr);
1507   }
1508   if (network->Jv) {
1509     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
1510     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
1511   }
1512 
1513   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
1514   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
1515   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
1516   if (network->vertex.sf) {
1517     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
1518   }
1519   /* edge */
1520   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
1521   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
1522   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
1523   if (network->edge.sf) {
1524     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
1525   }
1526   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
1527   network->edges = NULL;
1528   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
1529   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
1530 
1531   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
1532   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
1533   ierr = PetscFree(network->header);CHKERRQ(ierr);
1534   ierr = PetscFree(network);CHKERRQ(ierr);
1535   PetscFunctionReturn(0);
1536 }
1537 
1538 PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
1539 {
1540   PetscErrorCode ierr;
1541   DM_Network     *network = (DM_Network*) dm->data;
1542 
1543   PetscFunctionBegin;
1544   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
1545   PetscFunctionReturn(0);
1546 }
1547 
1548 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
1549 {
1550   PetscErrorCode ierr;
1551   DM_Network     *network = (DM_Network*) dm->data;
1552 
1553   PetscFunctionBegin;
1554   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
1555   PetscFunctionReturn(0);
1556 }
1557 
1558 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
1559 {
1560   PetscErrorCode ierr;
1561   DM_Network     *network = (DM_Network*) dm->data;
1562 
1563   PetscFunctionBegin;
1564   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
1565   PetscFunctionReturn(0);
1566 }
1567 
1568 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
1569 {
1570   PetscErrorCode ierr;
1571   DM_Network     *network = (DM_Network*) dm->data;
1572 
1573   PetscFunctionBegin;
1574   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
1575   PetscFunctionReturn(0);
1576 }
1577 
1578 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
1579 {
1580   PetscErrorCode ierr;
1581   DM_Network     *network = (DM_Network*) dm->data;
1582 
1583   PetscFunctionBegin;
1584   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
1585   PetscFunctionReturn(0);
1586 }
1587