xref: /petsc/src/dm/impls/network/network.c (revision 5cf7da580ae66613b1c0041b5e4fa11d82ae7164)
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   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);
271 
272   header->size[header->ndata] = component->size;
273   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
274   header->key[header->ndata] = componentkey;
275   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
276 
277   cvalue->data[header->ndata] = (void*)compvalue;
278   header->ndata++;
279   PetscFunctionReturn(0);
280 }
281 
282 #undef __FUNCT__
283 #define __FUNCT__ "DMNetworkGetNumComponents"
284 /*@
285   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
286 
287   Not Collective
288 
289   Input Parameters:
290 + dm - The DMNetwork object
291 . p  - vertex/edge point
292 
293   Output Parameters:
294 . numcomponents - Number of components at the vertex/edge
295 
296   Level: intermediate
297 
298 .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
299 @*/
300 PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
301 {
302   PetscErrorCode ierr;
303   PetscInt       offset;
304   DM_Network     *network = (DM_Network*)dm->data;
305 
306   PetscFunctionBegin;
307   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
308   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
309   PetscFunctionReturn(0);
310 }
311 
312 #undef __FUNCT__
313 #define __FUNCT__ "DMNetworkGetComponentTypeOffset"
314 /*@
315   DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the
316                                     component value from the component data array
317 
318   Not Collective
319 
320   Input Parameters:
321 + dm      - The DMNetwork object
322 . p       - vertex/edge point
323 - compnum - component number
324 
325   Output Parameters:
326 + compkey - the key obtained when registering the component
327 - offset  - offset into the component data array associated with the vertex/edge point
328 
329   Notes:
330   Typical usage:
331 
332   DMNetworkGetComponentDataArray(dm, &arr);
333   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
334   Loop over vertices or edges
335     DMNetworkGetNumComponents(dm,v,&numcomps);
336     Loop over numcomps
337       DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset);
338       compdata = (UserCompDataType)(arr+offset);
339 
340   Level: intermediate
341 
342 .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
343 @*/
344 PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
345 {
346   PetscErrorCode           ierr;
347   PetscInt                 offsetp;
348   DMNetworkComponentHeader header;
349   DM_Network               *network = (DM_Network*)dm->data;
350 
351   PetscFunctionBegin;
352   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
353   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
354   if (compkey) *compkey = header->key[compnum];
355   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
356   PetscFunctionReturn(0);
357 }
358 
359 #undef __FUNCT__
360 #define __FUNCT__ "DMNetworkGetVariableOffset"
361 /*@
362   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
363 
364   Not Collective
365 
366   Input Parameters:
367 + dm     - The DMNetwork object
368 - p      - the edge/vertex point
369 
370   Output Parameters:
371 . offset - the offset
372 
373   Level: intermediate
374 
375 .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
376 @*/
377 PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
378 {
379   PetscErrorCode ierr;
380   DM_Network     *network = (DM_Network*)dm->data;
381 
382   PetscFunctionBegin;
383   ierr = PetscSectionGetOffset(network->DofSection,p,offset);CHKERRQ(ierr);
384   PetscFunctionReturn(0);
385 }
386 
387 #undef __FUNCT__
388 #define __FUNCT__ "DMNetworkGetVariableGlobalOffset"
389 /*@
390   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
391 
392   Not Collective
393 
394   Input Parameters:
395 + dm      - The DMNetwork object
396 - p       - the edge/vertex point
397 
398   Output Parameters:
399 . offsetg - the offset
400 
401   Level: intermediate
402 
403 .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
404 @*/
405 PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
406 {
407   PetscErrorCode ierr;
408   DM_Network     *network = (DM_Network*)dm->data;
409 
410   PetscFunctionBegin;
411   ierr = PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);CHKERRQ(ierr);
412   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost node */
413   PetscFunctionReturn(0);
414 }
415 
416 #undef __FUNCT__
417 #define __FUNCT__ "DMNetworkAddNumVariables"
418 /*@
419   DMNetworkAddNumVariables - Add number of variables associated with a given point.
420 
421   Not Collective
422 
423   Input Parameters:
424 + dm   - The DMNetworkObject
425 . p    - the vertex/edge point
426 - nvar - number of additional variables
427 
428   Level: intermediate
429 
430 .seealso: DMNetworkSetNumVariables
431 @*/
432 PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
433 {
434   PetscErrorCode ierr;
435   DM_Network     *network = (DM_Network*)dm->data;
436 
437   PetscFunctionBegin;
438   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
439   PetscFunctionReturn(0);
440 }
441 
442 #undef __FUNCT__
443 #define __FUNCT__ "DMNetworkGetNumVariables"
444 /*@
445   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
446 
447   Not Collective
448 
449   Input Parameters:
450 + dm   - The DMNetworkObject
451 - p    - the vertex/edge point
452 
453   Output Parameters:
454 . nvar - number of variables
455 
456   Level: intermediate
457 
458 .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
459 @*/
460 PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
461 {
462   PetscErrorCode ierr;
463   DM_Network     *network = (DM_Network*)dm->data;
464 
465   PetscFunctionBegin;
466   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
467   PetscFunctionReturn(0);
468 }
469 
470 #undef __FUNCT__
471 #define __FUNCT__ "DMNetworkSetNumVariables"
472 /*@
473   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
474 
475   Not Collective
476 
477   Input Parameters:
478 + dm   - The DMNetworkObject
479 . p    - the vertex/edge point
480 - nvar - number of variables
481 
482   Level: intermediate
483 
484 .seealso: DMNetworkAddNumVariables
485 @*/
486 PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
487 {
488   PetscErrorCode ierr;
489   DM_Network     *network = (DM_Network*)dm->data;
490 
491   PetscFunctionBegin;
492   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
493   PetscFunctionReturn(0);
494 }
495 
496 /* Sets up the array that holds the data for all components and its associated section. This
497    function is called during DMSetUp() */
498 #undef __FUNCT__
499 #define __FUNCT__ "DMNetworkComponentSetUp"
500 PetscErrorCode DMNetworkComponentSetUp(DM dm)
501 {
502   PetscErrorCode              ierr;
503   DM_Network     *network = (DM_Network*)dm->data;
504   PetscInt                    arr_size;
505   PetscInt                    p,offset,offsetp;
506   DMNetworkComponentHeader header;
507   DMNetworkComponentValue  cvalue;
508   DMNetworkComponentGenericDataType      *componentdataarray;
509   PetscInt ncomp, i;
510 
511   PetscFunctionBegin;
512   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
513   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
514   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
515   componentdataarray = network->componentdataarray;
516   for (p = network->pStart; p < network->pEnd; p++) {
517     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
518     /* Copy header */
519     header = &network->header[p];
520     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
521     /* Copy data */
522     cvalue = &network->cvalue[p];
523     ncomp = header->ndata;
524     for (i = 0; i < ncomp; i++) {
525       offset = offsetp + network->dataheadersize + header->offset[i];
526       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
527     }
528   }
529   PetscFunctionReturn(0);
530 }
531 
532 /* Sets up the section for dofs. This routine is called during DMSetUp() */
533 #undef __FUNCT__
534 #define __FUNCT__ "DMNetworkVariablesSetUp"
535 PetscErrorCode DMNetworkVariablesSetUp(DM dm)
536 {
537   PetscErrorCode ierr;
538   DM_Network     *network = (DM_Network*)dm->data;
539 
540   PetscFunctionBegin;
541   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
542   PetscFunctionReturn(0);
543 }
544 
545 #undef __FUNCT__
546 #define __FUNCT__ "DMNetworkGetComponentDataArray"
547 /*@C
548   DMNetworkGetComponentDataArray - Returns the component data array
549 
550   Not Collective
551 
552   Input Parameters:
553 . dm - The DMNetwork Object
554 
555   Output Parameters:
556 . componentdataarray - array that holds data for all components
557 
558   Level: intermediate
559 
560 .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents
561 @*/
562 PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
563 {
564   DM_Network     *network = (DM_Network*)dm->data;
565 
566   PetscFunctionBegin;
567   *componentdataarray = network->componentdataarray;
568   PetscFunctionReturn(0);
569 }
570 
571 #undef __FUNCT__
572 #define __FUNCT__ "DMNetworkDistribute"
573 /*@
574   DMNetworkDistribute - Distributes the network and moves associated component data.
575 
576   Collective
577 
578   Input Parameter:
579 + oldDM - the original DMNetwork object
580 - overlap - The overlap of partitions, 0 is the default
581 
582   Output Parameter:
583 . distDM - the distributed DMNetwork object
584 
585   Notes:
586   This routine should be called only when using multiple processors.
587 
588   Distributes the network with <overlap>-overlapping partitioning of the edges.
589 
590   Level: intermediate
591 
592 .seealso: DMNetworkCreate
593 @*/
594 PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM)
595 {
596   PetscErrorCode ierr;
597   DM_Network     *oldDMnetwork = (DM_Network*)oldDM->data;
598   PetscSF        pointsf;
599   DM             newDM;
600   DM_Network     *newDMnetwork;
601 
602   PetscFunctionBegin;
603   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr);
604   newDMnetwork = (DM_Network*)newDM->data;
605   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
606   /* Distribute plex dm and dof section */
607   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
608   /* Distribute dof section */
609   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr);
610   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
611   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr);
612   /* Distribute data and associated section */
613   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
614   /* Destroy point SF */
615   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
616 
617   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
618   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
619   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
620   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
621   newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart;
622   newDMnetwork->NNodes = oldDMnetwork->NNodes;
623   newDMnetwork->NEdges = oldDMnetwork->NEdges;
624   /* Set Dof section as the default section for dm */
625   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
626   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
627 
628   *distDM = newDM;
629   PetscFunctionReturn(0);
630 }
631 
632 #undef __FUNCT__
633 #define __FUNCT__ "DMNetworkGetSupportingEdges"
634 /*@C
635   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
636 
637   Not Collective
638 
639   Input Parameters:
640 + dm - The DMNetwork object
641 - p  - the vertex point
642 
643   Output Paramters:
644 + nedges - number of edges connected to this vertex point
645 - edges  - List of edge points
646 
647   Level: intermediate
648 
649   Fortran Notes:
650   Since it returns an array, this routine is only available in Fortran 90, and you must
651   include petsc.h90 in your code.
652 
653 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
654 @*/
655 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
656 {
657   PetscErrorCode ierr;
658   DM_Network     *network = (DM_Network*)dm->data;
659 
660   PetscFunctionBegin;
661   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
662   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
663   PetscFunctionReturn(0);
664 }
665 
666 #undef __FUNCT__
667 #define __FUNCT__ "DMNetworkGetConnectedNodes"
668 /*@C
669   DMNetworkGetConnectedNodes - Return the connected vertices for this edge point
670 
671   Not Collective
672 
673   Input Parameters:
674 + dm - The DMNetwork object
675 - p  - the edge point
676 
677   Output Paramters:
678 . vertices  - vertices connected to this edge
679 
680   Level: intermediate
681 
682   Fortran Notes:
683   Since it returns an array, this routine is only available in Fortran 90, and you must
684   include petsc.h90 in your code.
685 
686 .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
687 @*/
688 PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
689 {
690   PetscErrorCode ierr;
691   DM_Network     *network = (DM_Network*)dm->data;
692 
693   PetscFunctionBegin;
694   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
695   PetscFunctionReturn(0);
696 }
697 
698 #undef __FUNCT__
699 #define __FUNCT__ "DMNetworkIsGhostVertex"
700 /*@
701   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
702 
703   Not Collective
704 
705   Input Parameters:
706 + dm - The DMNetwork object
707 . p  - the vertex point
708 
709   Output Parameter:
710 . isghost - TRUE if the vertex is a ghost point
711 
712   Level: intermediate
713 
714 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
715 @*/
716 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
717 {
718   PetscErrorCode ierr;
719   DM_Network     *network = (DM_Network*)dm->data;
720   PetscInt       offsetg;
721   PetscSection   sectiong;
722 
723   PetscFunctionBegin;
724   *isghost = PETSC_FALSE;
725   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
726   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
727   if (offsetg < 0) *isghost = PETSC_TRUE;
728   PetscFunctionReturn(0);
729 }
730 
731 #undef __FUNCT__
732 #define __FUNCT__ "DMSetUp_Network"
733 PetscErrorCode DMSetUp_Network(DM dm)
734 {
735   PetscErrorCode ierr;
736   DM_Network     *network=(DM_Network*)dm->data;
737 
738   PetscFunctionBegin;
739   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
740   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
741 
742   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
743   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
744   PetscFunctionReturn(0);
745 }
746 
747 #undef __FUNCT__
748 #define __FUNCT__ "DMNetworkHasJacobian"
749 /*@
750     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
751                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
752 
753     Collective
754 
755     Input Parameters:
756 +   dm - The DMNetwork object
757 .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
758 -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
759 
760     Level: intermediate
761 
762 @*/
763 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
764 {
765   DM_Network     *network=(DM_Network*)dm->data;
766 
767   PetscFunctionBegin;
768   network->userEdgeJacobian   = eflg;
769   network->userVertexJacobian = vflg;
770   PetscFunctionReturn(0);
771 }
772 
773 #undef __FUNCT__
774 #define __FUNCT__ "DMNetworkEdgeSetMatrix"
775 /*@
776     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
777 
778     Not Collective
779 
780     Input Parameters:
781 +   dm - The DMNetwork object
782 .   p  - the edge point
783 -   J - array (size = 3) of Jacobian submatrices for this edge point:
784         J[0]: this edge
785         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes()
786 
787     Level: intermediate
788 
789 .seealso: DMNetworkVertexSetMatrix
790 @*/
791 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
792 {
793   PetscErrorCode ierr;
794   DM_Network     *network=(DM_Network*)dm->data;
795 
796   PetscFunctionBegin;
797   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
798   if (!network->Je) {
799     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
800   }
801   network->Je[3*p]   = J[0];
802   network->Je[3*p+1] = J[1];
803   network->Je[3*p+2] = J[2];
804   PetscFunctionReturn(0);
805 }
806 
807 #undef __FUNCT__
808 #define __FUNCT__ "DMNetworkVertexSetMatrix"
809 /*@
810     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
811 
812     Not Collective
813 
814     Input Parameters:
815 +   dm - The DMNetwork object
816 .   p  - the vertex point
817 -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
818         J[0]:       this vertex
819         J[1+2*i]:   i-th supporting edge
820         J[1+2*i+1]: i-th connected vertex
821 
822     Level: intermediate
823 
824 .seealso: DMNetworkEdgeSetMatrix
825 @*/
826 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
827 {
828   PetscErrorCode ierr;
829   DM_Network     *network=(DM_Network*)dm->data;
830   PetscInt       i,*vptr,nedges,vStart,vEnd;
831   const PetscInt *edges;
832 
833   PetscFunctionBegin;
834   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
835 
836   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
837 
838   if (!network->Jv) {
839     PetscInt nNodes = network->nNodes,nedges_total;
840     /* count nvertex_total */
841     nedges_total = 0;
842     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
843     ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr);
844 
845     vptr[0] = 0;
846     for (i=0; i<nNodes; i++) {
847       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
848       nedges_total += nedges;
849       vptr[i+1] = vptr[i] + 2*nedges + 1;
850     }
851 
852     ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr);
853     network->Jvptr = vptr;
854   }
855 
856   vptr = network->Jvptr;
857   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
858 
859   /* Set Jacobian for each supporting edge and connected vertex */
860   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
861   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
862   PetscFunctionReturn(0);
863 }
864 
865 #undef __FUNCT__
866 #define __FUNCT__ "MatSetPreallocationDenseblock_private"
867 PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
868 {
869   PetscErrorCode ierr;
870   PetscInt       j;
871   PetscScalar    val=(PetscScalar)ncols;
872 
873   PetscFunctionBegin;
874   if (!ghost) {
875     for (j=0; j<nrows; j++) {
876       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
877     }
878   } else {
879     for (j=0; j<nrows; j++) {
880       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
881     }
882   }
883   PetscFunctionReturn(0);
884 }
885 
886 #undef __FUNCT__
887 #define __FUNCT__ "MatSetPreallocationUserblock_private"
888 PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
889 {
890   PetscErrorCode ierr;
891   PetscInt       j,ncols_u;
892   PetscScalar    val;
893 
894   PetscFunctionBegin;
895   if (!ghost) {
896     for (j=0; j<nrows; j++) {
897       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
898       val = (PetscScalar)ncols_u;
899       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
900       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
901     }
902   } else {
903     for (j=0; j<nrows; j++) {
904       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
905       val = (PetscScalar)ncols_u;
906       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
907       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
908     }
909   }
910   PetscFunctionReturn(0);
911 }
912 
913 #undef __FUNCT__
914 #define __FUNCT__ "MatSetPreallocationblock_private"
915 PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
916 {
917   PetscErrorCode ierr;
918 
919   PetscFunctionBegin;
920   if (Ju) {
921     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
922   } else {
923     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
924   }
925   PetscFunctionReturn(0);
926 }
927 
928 #undef __FUNCT__
929 #define __FUNCT__ "MatSetDenseblock_private"
930 PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
931 {
932   PetscErrorCode ierr;
933   PetscInt       j,*cols;
934   PetscScalar    *zeros;
935 
936   PetscFunctionBegin;
937   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
938   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
939   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
940   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
941   PetscFunctionReturn(0);
942 }
943 
944 #undef __FUNCT__
945 #define __FUNCT__ "MatSetUserblock_private"
946 PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
947 {
948   PetscErrorCode ierr;
949   PetscInt       j,M,N,row,col,ncols_u;
950   const PetscInt *cols;
951   PetscScalar    zero=0.0;
952 
953   PetscFunctionBegin;
954   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
955   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);
956 
957   for (row=0; row<nrows; row++) {
958     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
959     for (j=0; j<ncols_u; j++) {
960       col = cols[j] + cstart;
961       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
962     }
963     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
964   }
965   PetscFunctionReturn(0);
966 }
967 
968 #undef __FUNCT__
969 #define __FUNCT__ "MatSetblock_private"
970 PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
971 {
972   PetscErrorCode ierr;
973 
974   PetscFunctionBegin;
975   if (Ju) {
976     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
977   } else {
978     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
979   }
980   PetscFunctionReturn(0);
981 }
982 
983 #undef __FUNCT__
984 #define __FUNCT__ "DMCreateMatrix_Network"
985 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
986 {
987   PetscErrorCode ierr;
988   DM_Network     *network = (DM_Network*) dm->data;
989   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
990   PetscInt       cstart,ncols,j,e,v;
991   PetscBool      ghost,ghost_vc;
992   Mat            Juser;
993   PetscSection   sectionGlobal;
994   PetscInt       nedges,*vptr=NULL,vc;
995   const PetscInt *edges,*cone;
996   MPI_Comm       comm;
997   Vec            vd_nz,vo_nz;
998   PetscInt       *dnnz,*onnz;
999   PetscScalar    *vdnz,*vonz;
1000 
1001   PetscFunctionBegin;
1002   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1003     /* user does not provide Jacobian blocks */
1004     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1005     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1006     PetscFunctionReturn(0);
1007   }
1008 
1009   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
1010   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1011   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1012   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1013 
1014   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
1015   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
1016 
1017   /* (1) Set matrix preallocation */
1018   /*------------------------------*/
1019   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1020   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1021   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1022   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1023   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1024   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1025 
1026   /* Set preallocation for edges */
1027   /*-----------------------------*/
1028   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1029 
1030   for (e=eStart; e<eEnd; e++) {
1031     /* Get row indices */
1032     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1033     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1034     if (nrows) {
1035       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1036       ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
1037       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1038 
1039       /* Set preallocation for conntected vertices */
1040       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1041       for (v=0; v<2; v++) {
1042         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1043         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1044 
1045         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1046         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
1047         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1048       }
1049 
1050       /* Set preallocation for edge self */
1051       cstart = rstart;
1052       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1053       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1054       ierr = PetscFree(rows);CHKERRQ(ierr);
1055     }
1056   }
1057 
1058   /* Set preallocation for vertices */
1059   /*--------------------------------*/
1060   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1061   if (vEnd - vStart) {
1062     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1063     vptr = network->Jvptr;
1064     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1065   }
1066 
1067   for (v=vStart; v<vEnd; v++) {
1068     /* Get row indices */
1069     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1070     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1071     if (!nrows) continue;
1072 
1073     ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
1074     for (j=0; j<nrows; j++) rows[j] = j + rstart;
1075 
1076     /* Get supporting edges and connected vertices */
1077     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1078     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1079 
1080     for (e=0; e<nedges; e++) {
1081       /* Supporting edges */
1082       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1083       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1084 
1085       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1086       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1087 
1088       /* Connected vertices */
1089       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1090       vc = (v == cone[0]) ? cone[1]:cone[0];
1091       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1092 
1093       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
1094       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1095 
1096       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1097       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,(ghost_vc||ghost),vd_nz,vo_nz);CHKERRQ(ierr);
1098     }
1099 
1100     /* Set preallocation for vertex self */
1101     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1102     if (!ghost) {
1103       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1104       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1105       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1106     }
1107     ierr = PetscFree(rows);CHKERRQ(ierr);
1108   }
1109 
1110   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1111   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
1112 
1113   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
1114 
1115   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1116   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1117 
1118   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1119   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1120   for (j=0; j<localSize; j++) {
1121     dnnz[j] = (PetscInt)vdnz[j];
1122     onnz[j] = (PetscInt)vonz[j];
1123   }
1124   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1125   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1126   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1127   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1128 
1129   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
1130   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
1131   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1132 
1133   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
1134 
1135   /* (2) Set matrix entries for edges */
1136   /*----------------------------------*/
1137   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1138 
1139   for (e=eStart; e<eEnd; e++) {
1140     /* Get row indices */
1141     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1142     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1143     if (nrows) {
1144       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1145       ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
1146       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1147 
1148       /* Set matrix entries for conntected vertices */
1149       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1150       for (v=0; v<2; v++) {
1151         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1152         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1153 
1154         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1155         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1156       }
1157 
1158       /* Set matrix entries for edge self */
1159       cstart = rstart;
1160       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1161       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1162       ierr = PetscFree(rows);CHKERRQ(ierr);
1163     }
1164   }
1165 
1166   /* Set matrix entries for vertices */
1167   /*---------------------------------*/
1168   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1169 
1170     for (v=vStart; v<vEnd; v++) {
1171       /* Get row indices */
1172       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1173       ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1174       if (!nrows) continue;
1175 
1176       ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
1177       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1178 
1179       /* Get supporting edges and connected vertices */
1180       ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1181 
1182       for (e=0; e<nedges; e++) {
1183         /* Supporting edges */
1184         ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1185         ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1186 
1187         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1188         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1189 
1190         /* Connected vertices */
1191         ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1192         vc = (v == cone[0]) ? cone[1]:cone[0];
1193 
1194         ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
1195         ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1196 
1197         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1198       ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1199       }
1200 
1201       /* Set matrix entries for vertex self */
1202       ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1203       if (!ghost) {
1204         ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1205         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1206         ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1207       }
1208       ierr = PetscFree(rows);CHKERRQ(ierr);
1209     }
1210   }
1211   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1212   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1213 
1214   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1215   PetscFunctionReturn(0);
1216 }
1217 
1218 #undef __FUNCT__
1219 #define __FUNCT__ "DMDestroy_Network"
1220 PetscErrorCode DMDestroy_Network(DM dm)
1221 {
1222   PetscErrorCode ierr;
1223   DM_Network     *network = (DM_Network*) dm->data;
1224 
1225   PetscFunctionBegin;
1226   if (--network->refct > 0) PetscFunctionReturn(0);
1227   if (network->Je) {
1228     ierr = PetscFree(network->Je);CHKERRQ(ierr);
1229   }
1230   if (network->Jv) {
1231     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
1232     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
1233   }
1234   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
1235   network->edges = NULL;
1236   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
1237   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
1238 
1239   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
1240   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
1241   ierr = PetscFree(network->header);CHKERRQ(ierr);
1242   ierr = PetscFree(network);CHKERRQ(ierr);
1243   PetscFunctionReturn(0);
1244 }
1245 
1246 #undef __FUNCT__
1247 #define __FUNCT__ "DMView_Network"
1248 PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
1249 {
1250   PetscErrorCode ierr;
1251   DM_Network     *network = (DM_Network*) dm->data;
1252 
1253   PetscFunctionBegin;
1254   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
1255   PetscFunctionReturn(0);
1256 }
1257 
1258 #undef __FUNCT__
1259 #define __FUNCT__ "DMGlobalToLocalBegin_Network"
1260 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
1261 {
1262   PetscErrorCode ierr;
1263   DM_Network     *network = (DM_Network*) dm->data;
1264 
1265   PetscFunctionBegin;
1266   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
1267   PetscFunctionReturn(0);
1268 }
1269 
1270 #undef __FUNCT__
1271 #define __FUNCT__ "DMGlobalToLocalEnd_Network"
1272 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
1273 {
1274   PetscErrorCode ierr;
1275   DM_Network     *network = (DM_Network*) dm->data;
1276 
1277   PetscFunctionBegin;
1278   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
1279   PetscFunctionReturn(0);
1280 }
1281 
1282 #undef __FUNCT__
1283 #define __FUNCT__ "DMLocalToGlobalBegin_Network"
1284 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
1285 {
1286   PetscErrorCode ierr;
1287   DM_Network     *network = (DM_Network*) dm->data;
1288 
1289   PetscFunctionBegin;
1290   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
1291   PetscFunctionReturn(0);
1292 }
1293 
1294 #undef __FUNCT__
1295 #define __FUNCT__ "DMLocalToGlobalEnd_Network"
1296 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
1297 {
1298   PetscErrorCode ierr;
1299   DM_Network     *network = (DM_Network*) dm->data;
1300 
1301   PetscFunctionBegin;
1302   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
1303   PetscFunctionReturn(0);
1304 }
1305