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