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