xref: /petsc/src/dm/impls/network/network.c (revision ffa9b3b1b4df436c9af82c61ab815cba07b7b4bc)
1 #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
2 #include <petscdmplex.h>
3 #include <petscsf.h>
4 
5 #undef __FUNCT__
6 #define __FUNCT__ "DMNetworkSetSizes"
7 /*@
8   DMNetworkSetSizes - Sets the local and global vertices and edges.
9 
10   Collective on DM
11 
12   Input Parameters:
13 + dm - the dm object
14 . nV - number of local vertices
15 . nE - number of local edges
16 . NV - number of global vertices (or PETSC_DETERMINE)
17 - NE - number of global edges (or PETSC_DETERMINE)
18 
19    Notes
20    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
21 
22    You cannot change the sizes once they have been set
23 
24    Level: intermediate
25 
26 .seealso: DMNetworkCreate
27 @*/
28 PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE)
29 {
30   PetscErrorCode ierr;
31   DM_Network     *network = (DM_Network*) dm->data;
32   PetscInt       a[2],b[2];
33 
34   PetscFunctionBegin;
35   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
36   if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4);
37   if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5);
38   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);
39   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);
40   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);
41   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);
42   if (NE < 0 || NV < 0) {
43     a[0] = nV; a[1] = nE;
44     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
45     NV = b[0]; NE = b[1];
46   }
47   network->nNodes = nV;
48   network->NNodes = NV;
49   network->nEdges = nE;
50   network->NEdges = NE;
51   PetscFunctionReturn(0);
52 }
53 
54 #undef __FUNCT__
55 #define __FUNCT__ "DMNetworkSetEdgeList"
56 /*@
57   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
58 
59   Logically collective on DM
60 
61   Input Parameters:
62 . edges - list of edges
63 
64   Notes:
65   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
66   not be destroyed before the call to DMNetworkLayoutSetUp
67 
68   Level: intermediate
69 
70 .seealso: DMNetworkCreate, DMNetworkSetSizes
71 @*/
72 PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[])
73 {
74   DM_Network *network = (DM_Network*) dm->data;
75 
76   PetscFunctionBegin;
77   network->edges = edgelist;
78   PetscFunctionReturn(0);
79 }
80 
81 #undef __FUNCT__
82 #define __FUNCT__ "DMNetworkLayoutSetUp"
83 /*@
84   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
85 
86   Collective on DM
87 
88   Input Parameters
89 . DM - the dmnetwork object
90 
91   Notes:
92   This routine should be called after the network sizes and edgelists have been provided. It creates
93   the bare layout of the network and sets up the network to begin insertion of components.
94 
95   All the components should be registered before calling this routine.
96 
97   Level: intermediate
98 
99 .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
100 @*/
101 PetscErrorCode DMNetworkLayoutSetUp(DM dm)
102 {
103   PetscErrorCode ierr;
104   DM_Network     *network = (DM_Network*) dm->data;
105   PetscInt       dim = 1; /* One dimensional network */
106   PetscInt       numCorners=2;
107   PetscInt       spacedim=2;
108   double         *vertexcoords=NULL;
109   PetscInt       i;
110   PetscInt       ndata;
111 
112   PetscFunctionBegin;
113   if (network->nNodes) {
114     ierr = PetscCalloc1(numCorners*network->nNodes,&vertexcoords);CHKERRQ(ierr);
115   }
116   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
117   if (network->nNodes) {
118     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
119   }
120   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
121   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
122   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
123 
124   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
125   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
126   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
127   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
128 
129   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
130   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
131   for (i = network->pStart; i < network->pEnd; i++) {
132     network->header[i].ndata = 0;
133     ndata = network->header[i].ndata;
134     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
135     network->header[i].offset[ndata] = 0;
136   }
137   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
138   PetscFunctionReturn(0);
139 }
140 
141 #undef __FUNCT__
142 #define __FUNCT__ "DMNetworkRegisterComponent"
143 /*@
144   DMNetworkRegisterComponent - Registers the network component
145 
146   Logically collective on DM
147 
148   Input Parameters
149 + dm   - the network object
150 . name - the component name
151 - size - the storage size in bytes for this component data
152 
153    Output Parameters
154 .   key - an integer key that defines the component
155 
156    Notes
157    This routine should be called by all processors before calling DMNetworkLayoutSetup().
158 
159    Level: intermediate
160 
161 .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
162 @*/
163 PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
164 {
165   PetscErrorCode        ierr;
166   DM_Network            *network = (DM_Network*) dm->data;
167   DMNetworkComponent    *component=&network->component[network->ncomponent];
168   PetscBool             flg=PETSC_FALSE;
169   PetscInt              i;
170 
171   PetscFunctionBegin;
172 
173   for (i=0; i < network->ncomponent; i++) {
174     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
175     if (flg) {
176       *key = i;
177       PetscFunctionReturn(0);
178     }
179   }
180 
181   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
182   component->size = size/sizeof(DMNetworkComponentGenericDataType);
183   *key = network->ncomponent;
184   network->ncomponent++;
185   PetscFunctionReturn(0);
186 }
187 
188 #undef __FUNCT__
189 #define __FUNCT__ "DMNetworkGetVertexRange"
190 /*@
191   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
192 
193   Not Collective
194 
195   Input Parameters:
196 + dm - The DMNetwork object
197 
198   Output Paramters:
199 + vStart - The first vertex point
200 - vEnd   - One beyond the last vertex point
201 
202   Level: intermediate
203 
204 .seealso: DMNetworkGetEdgeRange
205 @*/
206 PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
207 {
208   DM_Network     *network = (DM_Network*)dm->data;
209 
210   PetscFunctionBegin;
211   if (vStart) *vStart = network->vStart;
212   if (vEnd) *vEnd = network->vEnd;
213   PetscFunctionReturn(0);
214 }
215 
216 #undef __FUNCT__
217 #define __FUNCT__ "DMNetworkGetEdgeRange"
218 /*@
219   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
220 
221   Not Collective
222 
223   Input Parameters:
224 + dm - The DMNetwork object
225 
226   Output Paramters:
227 + eStart - The first edge point
228 - eEnd   - One beyond the last edge point
229 
230   Level: intermediate
231 
232 .seealso: DMNetworkGetVertexRange
233 @*/
234 PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
235 {
236   DM_Network     *network = (DM_Network*)dm->data;
237 
238   PetscFunctionBegin;
239   if (eStart) *eStart = network->eStart;
240   if (eEnd) *eEnd = network->eEnd;
241   PetscFunctionReturn(0);
242 }
243 
244 #undef __FUNCT__
245 #define __FUNCT__ "DMNetworkAddComponent"
246 /*@
247   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
248 
249   Not Collective
250 
251   Input Parameters:
252 + dm           - The DMNetwork object
253 . p            - vertex/edge point
254 . componentkey - component key returned while registering the component
255 - compvalue    - pointer to the data structure for the component
256 
257   Level: intermediate
258 
259 .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
260 @*/
261 PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
262 {
263   DM_Network               *network = (DM_Network*)dm->data;
264   DMNetworkComponent       *component = &network->component[componentkey];
265   DMNetworkComponentHeader header = &network->header[p];
266   DMNetworkComponentValue  cvalue = &network->cvalue[p];
267   PetscErrorCode           ierr;
268 
269   PetscFunctionBegin;
270   header->size[header->ndata] = component->size;
271   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
272   header->key[header->ndata] = componentkey;
273   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
274 
275   cvalue->data[header->ndata] = (void*)compvalue;
276   header->ndata++;
277   PetscFunctionReturn(0);
278 }
279 
280 #undef __FUNCT__
281 #define __FUNCT__ "DMNetworkGetNumComponents"
282 /*@
283   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
284 
285   Not Collective
286 
287   Input Parameters:
288 + dm - The DMNetwork object
289 . p  - vertex/edge point
290 
291   Output Parameters:
292 . numcomponents - Number of components at the vertex/edge
293 
294   Level: intermediate
295 
296 .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
297 @*/
298 PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
299 {
300   PetscErrorCode ierr;
301   PetscInt       offset;
302   DM_Network     *network = (DM_Network*)dm->data;
303 
304   PetscFunctionBegin;
305   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
306   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
307   PetscFunctionReturn(0);
308 }
309 
310 #undef __FUNCT__
311 #define __FUNCT__ "DMNetworkGetComponentTypeOffset"
312 /*@
313   DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the
314                                     component value from the component data array
315 
316   Not Collective
317 
318   Input Parameters:
319 + dm      - The DMNetwork object
320 . p       - vertex/edge point
321 - compnum - component number
322 
323   Output Parameters:
324 + compkey - the key obtained when registering the component
325 - offset  - offset into the component data array associated with the vertex/edge point
326 
327   Notes:
328   Typical usage:
329 
330   DMNetworkGetComponentDataArray(dm, &arr);
331   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
332   Loop over vertices or edges
333     DMNetworkGetNumComponents(dm,v,&numcomps);
334     Loop over numcomps
335       DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset);
336       compdata = (UserCompDataType)(arr+offset);
337 
338   Level: intermediate
339 
340 .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
341 @*/
342 PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
343 {
344   PetscErrorCode           ierr;
345   PetscInt                 offsetp;
346   DMNetworkComponentHeader header;
347   DM_Network               *network = (DM_Network*)dm->data;
348 
349   PetscFunctionBegin;
350   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
351   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
352   if (compkey) *compkey = header->key[compnum];
353   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
354   PetscFunctionReturn(0);
355 }
356 
357 #undef __FUNCT__
358 #define __FUNCT__ "DMNetworkGetVariableOffset"
359 /*@
360   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
361 
362   Not Collective
363 
364   Input Parameters:
365 + dm     - The DMNetwork object
366 - p      - the edge/vertex point
367 
368   Output Parameters:
369 . offset - the offset
370 
371   Level: intermediate
372 
373 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
374 @*/
375 PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
376 {
377   PetscErrorCode ierr;
378   DM_Network     *network = (DM_Network*)dm->data;
379 
380   PetscFunctionBegin;
381   ierr = PetscSectionGetOffset(network->DofSection,p,offset);CHKERRQ(ierr);
382   PetscFunctionReturn(0);
383 }
384 
385 #undef __FUNCT__
386 #define __FUNCT__ "DMNetworkGetVariableGlobalOffset"
387 /*@
388   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
389 
390   Not Collective
391 
392   Input Parameters:
393 + dm      - The DMNetwork object
394 - p       - the edge/vertex point
395 
396   Output Parameters:
397 . offsetg - the offset
398 
399   Level: intermediate
400 
401 .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
402 @*/
403 PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
404 {
405   PetscErrorCode ierr;
406   DM_Network     *network = (DM_Network*)dm->data;
407 
408   PetscFunctionBegin;
409   ierr = PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);CHKERRQ(ierr);
410   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost node */
411   PetscFunctionReturn(0);
412 }
413 
414 #undef __FUNCT__
415 #define __FUNCT__ "DMNetworkAddNumVariables"
416 /*@
417   DMNetworkAddNumVariables - Add number of variables associated with a given point.
418 
419   Not Collective
420 
421   Input Parameters:
422 + dm   - The DMNetworkObject
423 . p    - the vertex/edge point
424 - nvar - number of additional variables
425 
426   Level: intermediate
427 
428 .seealso: DMNetworkSetNumVariables
429 @*/
430 PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
431 {
432   PetscErrorCode ierr;
433   DM_Network     *network = (DM_Network*)dm->data;
434 
435   PetscFunctionBegin;
436   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
437   PetscFunctionReturn(0);
438 }
439 
440 #undef __FUNCT__
441 #define __FUNCT__ "DMNetworkGetNumVariables"
442 /*@
443   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
444 
445   Not Collective
446 
447   Input Parameters:
448 + dm   - The DMNetworkObject
449 - p    - the vertex/edge point
450 
451   Output Parameters:
452 . nvar - number of variables
453 
454   Level: intermediate
455 
456 .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
457 @*/
458 PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
459 {
460   PetscErrorCode ierr;
461   DM_Network     *network = (DM_Network*)dm->data;
462 
463   PetscFunctionBegin;
464   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
465   PetscFunctionReturn(0);
466 }
467 
468 #undef __FUNCT__
469 #define __FUNCT__ "DMNetworkSetNumVariables"
470 /*@
471   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
472 
473   Not Collective
474 
475   Input Parameters:
476 + dm   - The DMNetworkObject
477 . p    - the vertex/edge point
478 - nvar - number of variables
479 
480   Level: intermediate
481 
482 .seealso: DMNetworkAddNumVariables
483 @*/
484 PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
485 {
486   PetscErrorCode ierr;
487   DM_Network     *network = (DM_Network*)dm->data;
488 
489   PetscFunctionBegin;
490   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
491   PetscFunctionReturn(0);
492 }
493 
494 /* Sets up the array that holds the data for all components and its associated section. This
495    function is called during DMSetUp() */
496 #undef __FUNCT__
497 #define __FUNCT__ "DMNetworkComponentSetUp"
498 PetscErrorCode DMNetworkComponentSetUp(DM dm)
499 {
500   PetscErrorCode              ierr;
501   DM_Network     *network = (DM_Network*)dm->data;
502   PetscInt                    arr_size;
503   PetscInt                    p,offset,offsetp;
504   DMNetworkComponentHeader header;
505   DMNetworkComponentValue  cvalue;
506   DMNetworkComponentGenericDataType      *componentdataarray;
507   PetscInt ncomp, i;
508 
509   PetscFunctionBegin;
510   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
511   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
512   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
513   componentdataarray = network->componentdataarray;
514   for (p = network->pStart; p < network->pEnd; p++) {
515     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
516     /* Copy header */
517     header = &network->header[p];
518     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
519     /* Copy data */
520     cvalue = &network->cvalue[p];
521     ncomp = header->ndata;
522     for (i = 0; i < ncomp; i++) {
523       offset = offsetp + network->dataheadersize + header->offset[i];
524       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
525     }
526   }
527   PetscFunctionReturn(0);
528 }
529 
530 /* Sets up the section for dofs. This routine is called during DMSetUp() */
531 #undef __FUNCT__
532 #define __FUNCT__ "DMNetworkVariablesSetUp"
533 PetscErrorCode DMNetworkVariablesSetUp(DM dm)
534 {
535   PetscErrorCode ierr;
536   DM_Network     *network = (DM_Network*)dm->data;
537 
538   PetscFunctionBegin;
539   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
540   PetscFunctionReturn(0);
541 }
542 
543 #undef __FUNCT__
544 #define __FUNCT__ "DMNetworkGetComponentDataArray"
545 /*@C
546   DMNetworkGetComponentDataArray - Returns the component data array
547 
548   Not Collective
549 
550   Input Parameters:
551 . dm - The DMNetwork Object
552 
553   Output Parameters:
554 . componentdataarray - array that holds data for all components
555 
556   Level: intermediate
557 
558 .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents
559 @*/
560 PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
561 {
562   DM_Network     *network = (DM_Network*)dm->data;
563 
564   PetscFunctionBegin;
565   *componentdataarray = network->componentdataarray;
566   PetscFunctionReturn(0);
567 }
568 
569 #undef __FUNCT__
570 #define __FUNCT__ "DMNetworkDistribute"
571 /*@
572   DMNetworkDistribute - Distributes the network and moves associated component data.
573 
574   Collective
575 
576   Input Parameter:
577 + oldDM - the original DMNetwork object
578 - overlap - The overlap of partitions, 0 is the default
579 
580   Output Parameter:
581 . distDM - the distributed DMNetwork object
582 
583   Notes:
584   This routine should be called only when using multiple processors.
585 
586   Distributes the network with <overlap>-overlapping partitioning of the edges.
587 
588   Level: intermediate
589 
590 .seealso: DMNetworkCreate
591 @*/
592 PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM)
593 {
594   PetscErrorCode ierr;
595   DM_Network     *oldDMnetwork = (DM_Network*)oldDM->data;
596   PetscSF        pointsf;
597   DM             newDM;
598   DM_Network     *newDMnetwork;
599 
600   PetscFunctionBegin;
601   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr);
602   newDMnetwork = (DM_Network*)newDM->data;
603   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
604   /* Distribute plex dm and dof section */
605   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
606   /* Distribute dof section */
607   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr);
608   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
609   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr);
610   /* Distribute data and associated section */
611   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
612   /* Destroy point SF */
613   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
614 
615   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
616   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
617   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
618   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
619   newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart;
620   newDMnetwork->NNodes = oldDMnetwork->NNodes;
621   newDMnetwork->NEdges = oldDMnetwork->NEdges;
622   /* Set Dof section as the default section for dm */
623   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
624   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
625 
626   *distDM = newDM;
627   PetscFunctionReturn(0);
628 }
629 
630 #undef __FUNCT__
631 #define __FUNCT__ "DMNetworkGetSupportingEdges"
632 /*@C
633   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
634 
635   Not Collective
636 
637   Input Parameters:
638 + dm - The DMNetwork object
639 - p  - the vertex point
640 
641   Output Paramters:
642 + nedges - number of edges connected to this vertex point
643 - edges  - List of edge points
644 
645   Level: intermediate
646 
647   Fortran Notes:
648   Since it returns an array, this routine is only available in Fortran 90, and you must
649   include petsc.h90 in your code.
650 
651 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
652 @*/
653 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
654 {
655   PetscErrorCode ierr;
656   DM_Network     *network = (DM_Network*)dm->data;
657 
658   PetscFunctionBegin;
659   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
660   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
661   PetscFunctionReturn(0);
662 }
663 
664 #undef __FUNCT__
665 #define __FUNCT__ "DMNetworkGetConnectedNodes"
666 /*@C
667   DMNetworkGetConnectedNodes - Return the connected vertices for this edge point
668 
669   Not Collective
670 
671   Input Parameters:
672 + dm - The DMNetwork object
673 - p  - the edge point
674 
675   Output Paramters:
676 . vertices  - vertices connected to this edge
677 
678   Level: intermediate
679 
680   Fortran Notes:
681   Since it returns an array, this routine is only available in Fortran 90, and you must
682   include petsc.h90 in your code.
683 
684 .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
685 @*/
686 PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
687 {
688   PetscErrorCode ierr;
689   DM_Network     *network = (DM_Network*)dm->data;
690 
691   PetscFunctionBegin;
692   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
693   PetscFunctionReturn(0);
694 }
695 
696 #undef __FUNCT__
697 #define __FUNCT__ "DMNetworkIsGhostVertex"
698 /*@
699   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
700 
701   Not Collective
702 
703   Input Parameters:
704 + dm - The DMNetwork object
705 . p  - the vertex point
706 
707   Output Parameter:
708 . isghost - TRUE if the vertex is a ghost point
709 
710   Level: intermediate
711 
712 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
713 @*/
714 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
715 {
716   PetscErrorCode ierr;
717   DM_Network     *network = (DM_Network*)dm->data;
718   PetscInt       offsetg;
719   PetscSection   sectiong;
720 
721   PetscFunctionBegin;
722   *isghost = PETSC_FALSE;
723   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
724   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
725   if (offsetg < 0) *isghost = PETSC_TRUE;
726   PetscFunctionReturn(0);
727 }
728 
729 #undef __FUNCT__
730 #define __FUNCT__ "DMSetUp_Network"
731 PetscErrorCode DMSetUp_Network(DM dm)
732 {
733   PetscErrorCode ierr;
734   DM_Network     *network=(DM_Network*)dm->data;
735 
736   PetscFunctionBegin;
737   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
738   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
739 
740   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
741   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
742   PetscFunctionReturn(0);
743 }
744 
745 #undef __FUNCT__
746 #define __FUNCT__ "DMNetworkHasJacobian"
747 /*@
748     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
749                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
750 
751     Collective
752 
753     Input Parameters:
754 +   dm - The DMNetwork object
755 .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
756 -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
757 
758     Level: intermediate
759 
760 @*/
761 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
762 {
763   DM_Network     *network=(DM_Network*)dm->data;
764 
765   PetscFunctionBegin;
766   network->userEdgeJacobian   = eflg;
767   network->userVertexJacobian = vflg;
768   PetscFunctionReturn(0);
769 }
770 
771 #undef __FUNCT__
772 #define __FUNCT__ "DMNetworkEdgeSetMatrix"
773 /*@
774     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
775 
776     Not Collective
777 
778     Input Parameters:
779 +   dm - The DMNetwork object
780 .   p  - the edge point
781 -   J - array (size = 3) of Jacobian submatrices for this edge point:
782         J[0]: this edge
783         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes()
784 
785     Level: intermediate
786 
787 .seealso: DMNetworkVertexSetMatrix
788 @*/
789 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
790 {
791   PetscErrorCode ierr;
792   DM_Network     *network=(DM_Network*)dm->data;
793 
794   PetscFunctionBegin;
795   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
796   if (!network->Je) {
797     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
798   }
799   network->Je[3*p]   = J[0];
800   network->Je[3*p+1] = J[1];
801   network->Je[3*p+2] = J[2];
802   PetscFunctionReturn(0);
803 }
804 
805 #undef __FUNCT__
806 #define __FUNCT__ "DMNetworkVertexSetMatrix"
807 /*@
808     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
809 
810     Not Collective
811 
812     Input Parameters:
813 +   dm - The DMNetwork object
814 .   p  - the vertex point
815 -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
816         J[0]:       this vertex
817         J[1+2*i]:   i-th supporting edge
818         J[1+2*i+1]: i-th connected vertex
819 
820     Level: intermediate
821 
822 .seealso: DMNetworkEdgeSetMatrix
823 @*/
824 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
825 {
826   PetscErrorCode ierr;
827   DM_Network     *network=(DM_Network*)dm->data;
828   PetscInt       i,*vptr,nedges,vStart,vEnd;
829   const PetscInt *edges;
830 
831   PetscFunctionBegin;
832   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
833 
834   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
835 
836   if (!network->Jv) {
837     PetscInt nNodes = network->nNodes,nedges_total;
838     /* count nvertex_total */
839     nedges_total = 0;
840     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
841     ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr);
842 
843     vptr[0] = 0;
844     for (i=0; i<nNodes; i++) {
845       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
846       nedges_total += nedges;
847       vptr[i+1] = vptr[i] + 2*nedges + 1;
848     }
849 
850     ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr);
851     network->Jvptr = vptr;
852   }
853 
854   vptr = network->Jvptr;
855   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
856 
857   /* Set Jacobian for each supporting edge and connected vertex */
858   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
859   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
860   PetscFunctionReturn(0);
861 }
862 
863 #undef __FUNCT__
864 #define __FUNCT__ "MatSetDenseblock_private"
865 PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
866 {
867   PetscErrorCode ierr;
868   PetscInt       j,*cols;
869   PetscScalar    *zeros;
870 
871   PetscFunctionBegin;
872   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
873   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
874   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
875   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
876   PetscFunctionReturn(0);
877 }
878 
879 #undef __FUNCT__
880 #define __FUNCT__ "MatSetUserblock_private"
881 PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
882 {
883   PetscErrorCode ierr;
884   PetscInt       j,M,N,row,col,ncols_u;
885   const PetscInt *cols;
886   PetscScalar    zero=0.0;
887 
888   PetscFunctionBegin;
889   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
890   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);
891 
892   for (row=0; row<nrows; row++) {
893     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
894     for (j=0; j<ncols_u; j++) {
895       col = cols[j] + cstart;
896       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
897     }
898     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
899   }
900   PetscFunctionReturn(0);
901 }
902 
903 #undef __FUNCT__
904 #define __FUNCT__ "MatSetblock_private"
905 PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
906 {
907   PetscErrorCode ierr;
908 
909   PetscFunctionBegin;
910   if (Ju) {
911     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
912   } else {
913     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
914   }
915   PetscFunctionReturn(0);
916 }
917 
918 #undef __FUNCT__
919 #define __FUNCT__ "DMCreateMatrix_Network"
920 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
921 {
922   PetscErrorCode ierr;
923   DM_Network     *network = (DM_Network*) dm->data;
924   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
925   PetscInt       cstart,ncols,j,e,v,*dnz,*onz,*dnzu,*onzu;
926   PetscBool      ghost;
927   Mat            Juser;
928   PetscSection   sectionGlobal;
929   PetscInt       nedges,*vptr,vc;
930   const PetscInt *edges,*cone;
931 
932   PetscFunctionBegin;
933   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
934     /* user does not provide Jacobian blocks */
935     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
936     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
937     PetscFunctionReturn(0);
938   }
939 
940   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
941   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
942   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
943   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
944 
945   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
946   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
947 
948   /* Preallocation - submatrix for an element (edge/vertex) is allocated as a dense block,
949    see DMCreateMatrix_Plex() */
950   ierr = PetscCalloc4(localSize,&dnz,localSize,&onz,localSize,&dnzu,localSize,&onzu);CHKERRQ(ierr);
951   ierr = DMPlexPreallocateOperator(network->plex,1,dnz,onz,dnzu,onzu,*J,PETSC_FALSE);CHKERRQ(ierr);
952   ierr = PetscFree4(dnz,onz,dnzu,onzu);CHKERRQ(ierr);
953 
954   /* Set matrix entries for edges */
955   /*------------------------------*/
956   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
957 
958   for (e=eStart; e<eEnd; e++) {
959     /* Get row indices */
960     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
961     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
962     if (nrows) {
963       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
964       ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
965       for (j=0; j<nrows; j++) rows[j] = j + rstart;
966 
967       /* Set matrix entries for conntected vertices */
968       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
969       for (v=0; v<2; v++) {
970         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
971         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
972 
973         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
974         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
975       }
976 
977       /* Set matrix entries for edge self */
978       cstart = rstart;
979       Juser = network->Je[3*e]; /* Jacobian(e,e) */
980       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
981       ierr = PetscFree(rows);CHKERRQ(ierr);
982     }
983   }
984 
985   /* Set matrix entries for vertices */
986   /*---------------------------------*/
987   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
988   if (vEnd - vStart) {
989     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"User must provide Jv");
990     if (!(vptr = network->Jvptr)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"User must provide vptr");
991   }
992 
993   for (v=vStart; v<vEnd; v++) {
994     /* Get row indices */
995     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
996     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
997     if (!nrows) continue;
998 
999     ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
1000     for (j=0; j<nrows; j++) rows[j] = j + rstart;
1001 
1002     /* Get supporting edges and connected vertices */
1003     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1004 
1005     for (e=0; e<nedges; e++) {
1006       /* Supporting edges */
1007       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1008       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1009 
1010       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1011       ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1012 
1013       /* Connected vertices */
1014       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1015       vc = (v == cone[0]) ? cone[1]:cone[0];
1016 
1017       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
1018       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1019 
1020       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1021       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1022     }
1023 
1024     /* Set matrix entries for vertex self */
1025     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1026     if (!ghost) {
1027       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1028       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1029       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1030     }
1031     ierr = PetscFree(rows);CHKERRQ(ierr);
1032   }
1033   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1034   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1035 
1036   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1037   PetscFunctionReturn(0);
1038 }
1039 
1040 #undef __FUNCT__
1041 #define __FUNCT__ "DMDestroy_Network"
1042 PetscErrorCode DMDestroy_Network(DM dm)
1043 {
1044   PetscErrorCode ierr;
1045   DM_Network     *network = (DM_Network*) dm->data;
1046 
1047   PetscFunctionBegin;
1048   if (--network->refct > 0) PetscFunctionReturn(0);
1049   if (network->Je) {
1050     ierr = PetscFree(network->Je);CHKERRQ(ierr);
1051   }
1052   if (network->Jv) {
1053     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
1054     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
1055   }
1056   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
1057   network->edges = NULL;
1058   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
1059   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
1060 
1061   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
1062   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
1063   ierr = PetscFree(network->header);CHKERRQ(ierr);
1064   ierr = PetscFree(network);CHKERRQ(ierr);
1065   PetscFunctionReturn(0);
1066 }
1067 
1068 #undef __FUNCT__
1069 #define __FUNCT__ "DMView_Network"
1070 PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
1071 {
1072   PetscErrorCode ierr;
1073   DM_Network     *network = (DM_Network*) dm->data;
1074 
1075   PetscFunctionBegin;
1076   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
1077   PetscFunctionReturn(0);
1078 }
1079 
1080 #undef __FUNCT__
1081 #define __FUNCT__ "DMGlobalToLocalBegin_Network"
1082 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
1083 {
1084   PetscErrorCode ierr;
1085   DM_Network     *network = (DM_Network*) dm->data;
1086 
1087   PetscFunctionBegin;
1088   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
1089   PetscFunctionReturn(0);
1090 }
1091 
1092 #undef __FUNCT__
1093 #define __FUNCT__ "DMGlobalToLocalEnd_Network"
1094 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
1095 {
1096   PetscErrorCode ierr;
1097   DM_Network     *network = (DM_Network*) dm->data;
1098 
1099   PetscFunctionBegin;
1100   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
1101   PetscFunctionReturn(0);
1102 }
1103 
1104 #undef __FUNCT__
1105 #define __FUNCT__ "DMLocalToGlobalBegin_Network"
1106 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
1107 {
1108   PetscErrorCode ierr;
1109   DM_Network     *network = (DM_Network*) dm->data;
1110 
1111   PetscFunctionBegin;
1112   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
1113   PetscFunctionReturn(0);
1114 }
1115 
1116 #undef __FUNCT__
1117 #define __FUNCT__ "DMLocalToGlobalEnd_Network"
1118 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
1119 {
1120   PetscErrorCode ierr;
1121   DM_Network     *network = (DM_Network*) dm->data;
1122 
1123   PetscFunctionBegin;
1124   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
1125   PetscFunctionReturn(0);
1126 }
1127