xref: /petsc/src/dm/impls/network/network.c (revision 0731d606733b2d54a4589118097efd7d0c58699d)
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 
667   PetscFunctionReturn(0);
668 }
669 
670 #undef __FUNCT__
671 #define __FUNCT__ "DMNetworkAssembleGraphStructures"
672 /*@
673   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
674 
675   Collective
676 
677   Input Parameters:
678 . dm   - The DMNetworkObject
679 
680   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
681 
682   points = [0 1 2 3 4 5 6]
683 
684   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
685 
686   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
687 
688   Level: intermediate
689 
690 @*/
691 PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
692 {
693   PetscErrorCode ierr;
694   MPI_Comm               comm;
695   PetscMPIInt            rank, numProcs;
696   DM_Network     *network = (DM_Network*)dm->data;
697   PetscFunctionBegin;
698 
699   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
700   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
701   ierr = MPI_Comm_size(comm, &numProcs);CHKERRQ(ierr);
702 
703   /* Create maps for vertices and edges */
704   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
705   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
706 
707   /* Create local sub-sections */
708   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
709   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
710 
711   if (numProcs > 1) {
712     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
713     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
714   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
715   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
716   } else {
717   /* create structures for vertex */
718   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
719   /* create structures for edge */
720   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
721   }
722 
723 
724   /* Add viewers */
725   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
726   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
727   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
728   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
729 
730   PetscFunctionReturn(0);
731 }
732 #undef __FUNCT__
733 #define __FUNCT__ "DMNetworkDistribute"
734 /*@
735   DMNetworkDistribute - Distributes the network and moves associated component data.
736 
737   Collective
738 
739   Input Parameter:
740 + oldDM - the original DMNetwork object
741 - overlap - The overlap of partitions, 0 is the default
742 
743   Output Parameter:
744 . distDM - the distributed DMNetwork object
745 
746   Notes:
747   This routine should be called only when using multiple processors.
748 
749   Distributes the network with <overlap>-overlapping partitioning of the edges.
750 
751   Level: intermediate
752 
753 .seealso: DMNetworkCreate
754 @*/
755 PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM)
756 {
757   PetscErrorCode ierr;
758   DM_Network     *oldDMnetwork = (DM_Network*)oldDM->data;
759   PetscSF        pointsf;
760   DM             newDM;
761   DM_Network     *newDMnetwork;
762 
763   PetscFunctionBegin;
764   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr);
765   newDMnetwork = (DM_Network*)newDM->data;
766   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
767   /* Distribute plex dm and dof section */
768   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
769   /* Distribute dof section */
770   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr);
771   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
772   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr);
773   /* Distribute data and associated section */
774   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
775 
776 
777   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
778   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
779   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
780   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
781   newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart;
782   newDMnetwork->NNodes = oldDMnetwork->NNodes;
783   newDMnetwork->NEdges = oldDMnetwork->NEdges;
784 
785   /* Set Dof section as the default section for dm */
786   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
787   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
788 
789   /* Destroy point SF */
790   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
791 
792   *distDM = newDM;
793   PetscFunctionReturn(0);
794 }
795 
796 #undef __FUNCT__
797 #define __FUNCT__ "PetscSFGetSubSF"
798 /*@C
799   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
800 
801   Input Parameters:
802 + masterSF - the original SF structure
803 - map      - a ISLocalToGlobal mapping that contains the subset of points
804 
805   Output Parameters:
806 . subSF    - a subset of the masterSF for the desired subset.
807 */
808 
809 PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
810 
811   PetscErrorCode        ierr;
812   PetscInt              nroots, nleaves, *ilocal_sub;
813   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
814   PetscInt              *local_points, *remote_points;
815   PetscSFNode           *iremote_sub;
816   const PetscInt        *ilocal;
817   const PetscSFNode     *iremote;
818 
819   PetscFunctionBegin;
820   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
821 
822   /* Look for leaves that pertain to the subset of points. Get the local ordering */
823   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
824   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
825   for (i = 0; i < nleaves; i++) {
826     if (ilocal_map[i] != -1) nleaves_sub += 1;
827   }
828   /* Re-number ilocal with subset numbering. Need information from roots */
829   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
830   for (i = 0; i < nroots; i++) local_points[i] = i;
831   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
832   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
833   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
834   /* Fill up graph using local (that is, local to the subset) numbering. */
835   ierr = PetscMalloc2(nleaves_sub,&ilocal_sub,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] = iremote[i];
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_OWN_POINTER);CHKERRQ(ierr);
852   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
853 
854   PetscFunctionReturn(0);
855 }
856 
857 
858 #undef __FUNCT__
859 #define __FUNCT__ "DMNetworkGetSupportingEdges"
860 /*@C
861   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
862 
863   Not Collective
864 
865   Input Parameters:
866 + dm - The DMNetwork object
867 - p  - the vertex point
868 
869   Output Paramters:
870 + nedges - number of edges connected to this vertex point
871 - edges  - List of edge points
872 
873   Level: intermediate
874 
875   Fortran Notes:
876   Since it returns an array, this routine is only available in Fortran 90, and you must
877   include petsc.h90 in your code.
878 
879 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
880 @*/
881 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
882 {
883   PetscErrorCode ierr;
884   DM_Network     *network = (DM_Network*)dm->data;
885 
886   PetscFunctionBegin;
887   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
888   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
889   PetscFunctionReturn(0);
890 }
891 
892 #undef __FUNCT__
893 #define __FUNCT__ "DMNetworkGetConnectedNodes"
894 /*@C
895   DMNetworkGetConnectedNodes - Return the connected vertices for this edge point
896 
897   Not Collective
898 
899   Input Parameters:
900 + dm - The DMNetwork object
901 - p  - the edge point
902 
903   Output Paramters:
904 . vertices  - vertices connected to this edge
905 
906   Level: intermediate
907 
908   Fortran Notes:
909   Since it returns an array, this routine is only available in Fortran 90, and you must
910   include petsc.h90 in your code.
911 
912 .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
913 @*/
914 PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
915 {
916   PetscErrorCode ierr;
917   DM_Network     *network = (DM_Network*)dm->data;
918 
919   PetscFunctionBegin;
920   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
921   PetscFunctionReturn(0);
922 }
923 
924 #undef __FUNCT__
925 #define __FUNCT__ "DMNetworkIsGhostVertex"
926 /*@
927   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
928 
929   Not Collective
930 
931   Input Parameters:
932 + dm - The DMNetwork object
933 . p  - the vertex point
934 
935   Output Parameter:
936 . isghost - TRUE if the vertex is a ghost point
937 
938   Level: intermediate
939 
940 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
941 @*/
942 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
943 {
944   PetscErrorCode ierr;
945   DM_Network     *network = (DM_Network*)dm->data;
946   PetscInt       offsetg;
947   PetscSection   sectiong;
948 
949   PetscFunctionBegin;
950   *isghost = PETSC_FALSE;
951   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
952   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
953   if (offsetg < 0) *isghost = PETSC_TRUE;
954   PetscFunctionReturn(0);
955 }
956 
957 #undef __FUNCT__
958 #define __FUNCT__ "DMSetUp_Network"
959 PetscErrorCode DMSetUp_Network(DM dm)
960 {
961   PetscErrorCode ierr;
962   DM_Network     *network=(DM_Network*)dm->data;
963 
964   PetscFunctionBegin;
965   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
966   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
967 
968   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
969   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
970   PetscFunctionReturn(0);
971 }
972 
973 #undef __FUNCT__
974 #define __FUNCT__ "DMNetworkHasJacobian"
975 /*@
976     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
977                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
978 
979     Collective
980 
981     Input Parameters:
982 +   dm - The DMNetwork object
983 .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
984 -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
985 
986     Level: intermediate
987 
988 @*/
989 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
990 {
991   DM_Network     *network=(DM_Network*)dm->data;
992 
993   PetscFunctionBegin;
994   network->userEdgeJacobian   = eflg;
995   network->userVertexJacobian = vflg;
996   PetscFunctionReturn(0);
997 }
998 
999 #undef __FUNCT__
1000 #define __FUNCT__ "DMNetworkEdgeSetMatrix"
1001 /*@
1002     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
1003 
1004     Not Collective
1005 
1006     Input Parameters:
1007 +   dm - The DMNetwork object
1008 .   p  - the edge point
1009 -   J - array (size = 3) of Jacobian submatrices for this edge point:
1010         J[0]: this edge
1011         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes()
1012 
1013     Level: intermediate
1014 
1015 .seealso: DMNetworkVertexSetMatrix
1016 @*/
1017 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
1018 {
1019   PetscErrorCode ierr;
1020   DM_Network     *network=(DM_Network*)dm->data;
1021 
1022   PetscFunctionBegin;
1023   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
1024   if (!network->Je) {
1025     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
1026   }
1027   network->Je[3*p]   = J[0];
1028   network->Je[3*p+1] = J[1];
1029   network->Je[3*p+2] = J[2];
1030   PetscFunctionReturn(0);
1031 }
1032 
1033 #undef __FUNCT__
1034 #define __FUNCT__ "DMNetworkVertexSetMatrix"
1035 /*@
1036     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
1037 
1038     Not Collective
1039 
1040     Input Parameters:
1041 +   dm - The DMNetwork object
1042 .   p  - the vertex point
1043 -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
1044         J[0]:       this vertex
1045         J[1+2*i]:   i-th supporting edge
1046         J[1+2*i+1]: i-th connected vertex
1047 
1048     Level: intermediate
1049 
1050 .seealso: DMNetworkEdgeSetMatrix
1051 @*/
1052 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
1053 {
1054   PetscErrorCode ierr;
1055   DM_Network     *network=(DM_Network*)dm->data;
1056   PetscInt       i,*vptr,nedges,vStart,vEnd;
1057   const PetscInt *edges;
1058 
1059   PetscFunctionBegin;
1060   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1061 
1062   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1063 
1064   if (!network->Jv) {
1065     PetscInt nNodes = network->nNodes,nedges_total;
1066     /* count nvertex_total */
1067     nedges_total = 0;
1068     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1069     ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr);
1070 
1071     vptr[0] = 0;
1072     for (i=0; i<nNodes; i++) {
1073       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1074       nedges_total += nedges;
1075       vptr[i+1] = vptr[i] + 2*nedges + 1;
1076     }
1077 
1078     ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr);
1079     network->Jvptr = vptr;
1080   }
1081 
1082   vptr = network->Jvptr;
1083   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
1084 
1085   /* Set Jacobian for each supporting edge and connected vertex */
1086   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1087   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1088   PetscFunctionReturn(0);
1089 }
1090 
1091 #undef __FUNCT__
1092 #define __FUNCT__ "MatSetPreallocationDenseblock_private"
1093 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1094 {
1095   PetscErrorCode ierr;
1096   PetscInt       j;
1097   PetscScalar    val=(PetscScalar)ncols;
1098 
1099   PetscFunctionBegin;
1100   if (!ghost) {
1101     for (j=0; j<nrows; j++) {
1102       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1103     }
1104   } else {
1105     for (j=0; j<nrows; j++) {
1106       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1107     }
1108   }
1109   PetscFunctionReturn(0);
1110 }
1111 
1112 #undef __FUNCT__
1113 #define __FUNCT__ "MatSetPreallocationUserblock_private"
1114 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1115 {
1116   PetscErrorCode ierr;
1117   PetscInt       j,ncols_u;
1118   PetscScalar    val;
1119 
1120   PetscFunctionBegin;
1121   if (!ghost) {
1122     for (j=0; j<nrows; j++) {
1123       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1124       val = (PetscScalar)ncols_u;
1125       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1126       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1127     }
1128   } else {
1129     for (j=0; j<nrows; j++) {
1130       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1131       val = (PetscScalar)ncols_u;
1132       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1133       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1134     }
1135   }
1136   PetscFunctionReturn(0);
1137 }
1138 
1139 #undef __FUNCT__
1140 #define __FUNCT__ "MatSetPreallocationblock_private"
1141 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1142 {
1143   PetscErrorCode ierr;
1144 
1145   PetscFunctionBegin;
1146   if (Ju) {
1147     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
1148   } else {
1149     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
1150   }
1151   PetscFunctionReturn(0);
1152 }
1153 
1154 #undef __FUNCT__
1155 #define __FUNCT__ "MatSetDenseblock_private"
1156 PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1157 {
1158   PetscErrorCode ierr;
1159   PetscInt       j,*cols;
1160   PetscScalar    *zeros;
1161 
1162   PetscFunctionBegin;
1163   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1164   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1165   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1166   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
1167   PetscFunctionReturn(0);
1168 }
1169 
1170 #undef __FUNCT__
1171 #define __FUNCT__ "MatSetUserblock_private"
1172 PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1173 {
1174   PetscErrorCode ierr;
1175   PetscInt       j,M,N,row,col,ncols_u;
1176   const PetscInt *cols;
1177   PetscScalar    zero=0.0;
1178 
1179   PetscFunctionBegin;
1180   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
1181   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);
1182 
1183   for (row=0; row<nrows; row++) {
1184     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
1185     for (j=0; j<ncols_u; j++) {
1186       col = cols[j] + cstart;
1187       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
1188     }
1189     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
1190   }
1191   PetscFunctionReturn(0);
1192 }
1193 
1194 #undef __FUNCT__
1195 #define __FUNCT__ "MatSetblock_private"
1196 PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1197 {
1198   PetscErrorCode ierr;
1199 
1200   PetscFunctionBegin;
1201   if (Ju) {
1202     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1203   } else {
1204     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1205   }
1206   PetscFunctionReturn(0);
1207 }
1208 
1209 
1210 #undef __FUNCT__
1211 #define __FUNCT__ "CreateSubGlobalToLocalMapping_private"
1212 /* 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.
1213 */
1214 PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
1215 {
1216   PetscErrorCode ierr;
1217   PetscInt       i, size, dof;
1218   PetscInt       *glob2loc;
1219 
1220   PetscFunctionBegin;
1221   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
1222   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
1223 
1224   for (i = 0; i < size; i++) {
1225     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
1226     dof = (dof >= 0) ? dof : -(dof + 1);
1227     glob2loc[i] = dof;
1228   }
1229 
1230   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
1231 #if 0
1232   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1233 #endif
1234   PetscFunctionReturn(0);
1235 }
1236 
1237 #undef __FUNCT__
1238 #define __FUNCT__ "DMCreateMatrix_Network"
1239 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
1240 {
1241   PetscErrorCode ierr;
1242   PetscMPIInt    rank, size;
1243   DM_Network     *network = (DM_Network*) dm->data;
1244   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1245   PetscInt       cstart,ncols,j,e,v;
1246   PetscBool      ghost,ghost_vc,ghost2,isNest;
1247   Mat            Juser;
1248   PetscSection   sectionGlobal;
1249   PetscInt       nedges,*vptr=NULL,vc,*rows_v;
1250   const PetscInt *edges,*cone;
1251   MPI_Comm       comm;
1252   MatType        mtype;
1253   Vec            vd_nz,vo_nz;
1254   PetscInt       *dnnz,*onnz;
1255   PetscScalar    *vdnz,*vonz;
1256 
1257   PetscFunctionBegin;
1258   mtype = dm->mattype;
1259   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
1260 
1261   if (isNest) {
1262     /* ierr = DMCreateMatrix_Network_Nest(); */
1263     PetscInt   eDof, vDof;
1264     Mat        j11, j12, j21, j22, bA[2][2];
1265     ISLocalToGlobalMapping eISMap, vISMap;
1266 
1267     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1268     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1269     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1270 
1271     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
1272     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
1273 
1274 #if 0
1275     printf("rank[%d] vdof = %d. edof = %d\n", rank, vDof, eDof);
1276 #endif
1277 
1278     ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr);
1279     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1280     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
1281 
1282     ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr);
1283     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
1284     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
1285 
1286     ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr);
1287     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1288     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
1289 
1290     ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr);
1291     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1292     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
1293 
1294     bA[0][0] = j11;
1295     bA[0][1] = j12;
1296     bA[1][0] = j21;
1297     bA[1][1] = j22;
1298 
1299     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
1300     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
1301 
1302     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
1303     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
1304     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
1305     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
1306 
1307     ierr = MatSetUp(j11);CHKERRQ(ierr);
1308     ierr = MatSetUp(j12);CHKERRQ(ierr);
1309     ierr = MatSetUp(j21);CHKERRQ(ierr);
1310     ierr = MatSetUp(j22);CHKERRQ(ierr);
1311 
1312     ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
1313     ierr = MatSetUp(*J);CHKERRQ(ierr);
1314     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
1315     ierr = MatDestroy(&j11);CHKERRQ(ierr);
1316     ierr = MatDestroy(&j12);CHKERRQ(ierr);
1317     ierr = MatDestroy(&j21);CHKERRQ(ierr);
1318     ierr = MatDestroy(&j22);CHKERRQ(ierr);
1319 
1320     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1321     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1322 
1323     /* Free structures */
1324     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
1325     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
1326 
1327     PetscFunctionReturn(0);
1328   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1329     /* user does not provide Jacobian blocks */
1330     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1331     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1332     PetscFunctionReturn(0);
1333   }
1334 
1335   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
1336   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1337   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1338   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1339 
1340   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
1341   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
1342 
1343   /* (1) Set matrix preallocation */
1344   /*------------------------------*/
1345   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1346   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1347   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1348   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1349   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1350   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1351 
1352   /* Set preallocation for edges */
1353   /*-----------------------------*/
1354   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1355 
1356   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1357   for (e=eStart; e<eEnd; e++) {
1358     /* Get row indices */
1359     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1360     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1361     if (nrows) {
1362       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1363       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1364 
1365       /* Set preallocation for conntected vertices */
1366       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1367       for (v=0; v<2; v++) {
1368         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1369 
1370         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1371         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
1372         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1373       }
1374 
1375       /* Set preallocation for edge self */
1376       cstart = rstart;
1377       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1378       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1379     }
1380   }
1381 
1382   /* Set preallocation for vertices */
1383   /*--------------------------------*/
1384   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1385   if (vEnd - vStart) {
1386     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1387     vptr = network->Jvptr;
1388     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1389   }
1390 
1391   for (v=vStart; v<vEnd; v++) {
1392     /* Get row indices */
1393     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1394     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1395     if (!nrows) continue;
1396 
1397     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1398     if (ghost) {
1399       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1400     } else {
1401       rows_v = rows;
1402     }
1403 
1404     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1405 
1406     /* Get supporting edges and connected vertices */
1407     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1408 
1409     for (e=0; e<nedges; e++) {
1410       /* Supporting edges */
1411       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1412       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1413 
1414       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1415       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1416 
1417       /* Connected vertices */
1418       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1419       vc = (v == cone[0]) ? cone[1]:cone[0];
1420       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1421 
1422       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1423 
1424       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1425       if (ghost_vc||ghost) {
1426         ghost2 = PETSC_TRUE;
1427       } else {
1428         ghost2 = PETSC_FALSE;
1429       }
1430       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1431     }
1432 
1433     /* Set preallocation for vertex self */
1434     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1435     if (!ghost) {
1436       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1437       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1438       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1439     }
1440     if (ghost) {
1441       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1442     }
1443   }
1444 
1445   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1446   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
1447 
1448   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
1449 
1450   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1451   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1452 
1453   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1454   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1455   for (j=0; j<localSize; j++) {
1456     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1457     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1458   }
1459   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1460   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1461   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1462   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1463 
1464   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
1465   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
1466   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1467 
1468   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
1469 
1470   /* (2) Set matrix entries for edges */
1471   /*----------------------------------*/
1472   for (e=eStart; e<eEnd; e++) {
1473     /* Get row indices */
1474     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1475     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1476     if (nrows) {
1477       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1478 
1479       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1480 
1481       /* Set matrix entries for conntected vertices */
1482       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1483       for (v=0; v<2; v++) {
1484         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1485         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1486 
1487         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1488         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1489       }
1490 
1491       /* Set matrix entries for edge self */
1492       cstart = rstart;
1493       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1494       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1495     }
1496   }
1497 
1498   /* Set matrix entries for vertices */
1499   /*---------------------------------*/
1500   for (v=vStart; v<vEnd; v++) {
1501     /* Get row indices */
1502     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1503     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1504     if (!nrows) continue;
1505 
1506     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1507     if (ghost) {
1508       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1509     } else {
1510       rows_v = rows;
1511     }
1512     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1513 
1514     /* Get supporting edges and connected vertices */
1515     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1516 
1517     for (e=0; e<nedges; e++) {
1518       /* Supporting edges */
1519       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1520       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1521 
1522       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1523       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1524 
1525       /* Connected vertices */
1526       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1527       vc = (v == cone[0]) ? cone[1]:cone[0];
1528 
1529       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
1530       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1531 
1532       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1533       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1534     }
1535 
1536     /* Set matrix entries for vertex self */
1537     if (!ghost) {
1538       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1539       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1540       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1541     }
1542     if (ghost) {
1543       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1544     }
1545   }
1546   ierr = PetscFree(rows);CHKERRQ(ierr);
1547 
1548   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1549   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1550 
1551   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1552   PetscFunctionReturn(0);
1553 }
1554 
1555 #undef __FUNCT__
1556 #define __FUNCT__ "DMDestroy_Network"
1557 PetscErrorCode DMDestroy_Network(DM dm)
1558 {
1559   PetscErrorCode ierr;
1560   DM_Network     *network = (DM_Network*) dm->data;
1561 
1562   PetscFunctionBegin;
1563   /* vertex */
1564   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
1565   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
1566   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
1567   if (network->vertex.sf) {
1568     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
1569   }
1570   /* edge */
1571   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
1572   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
1573   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
1574   if (network->edge.sf) {
1575     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
1576   }
1577   if (--network->refct > 0) PetscFunctionReturn(0);
1578   if (network->Je) {
1579     ierr = PetscFree(network->Je);CHKERRQ(ierr);
1580   }
1581   if (network->Jv) {
1582     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
1583     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
1584   }
1585   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
1586   network->edges = NULL;
1587   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
1588   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
1589 
1590   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
1591   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
1592   ierr = PetscFree(network->header);CHKERRQ(ierr);
1593   ierr = PetscFree(network);CHKERRQ(ierr);
1594   PetscFunctionReturn(0);
1595 }
1596 
1597 #undef __FUNCT__
1598 #define __FUNCT__ "DMView_Network"
1599 PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
1600 {
1601   PetscErrorCode ierr;
1602   DM_Network     *network = (DM_Network*) dm->data;
1603 
1604   PetscFunctionBegin;
1605   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
1606   PetscFunctionReturn(0);
1607 }
1608 
1609 #undef __FUNCT__
1610 #define __FUNCT__ "DMGlobalToLocalBegin_Network"
1611 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
1612 {
1613   PetscErrorCode ierr;
1614   DM_Network     *network = (DM_Network*) dm->data;
1615 
1616   PetscFunctionBegin;
1617   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
1618   PetscFunctionReturn(0);
1619 }
1620 
1621 #undef __FUNCT__
1622 #define __FUNCT__ "DMGlobalToLocalEnd_Network"
1623 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
1624 {
1625   PetscErrorCode ierr;
1626   DM_Network     *network = (DM_Network*) dm->data;
1627 
1628   PetscFunctionBegin;
1629   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
1630   PetscFunctionReturn(0);
1631 }
1632 
1633 #undef __FUNCT__
1634 #define __FUNCT__ "DMLocalToGlobalBegin_Network"
1635 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
1636 {
1637   PetscErrorCode ierr;
1638   DM_Network     *network = (DM_Network*) dm->data;
1639 
1640   PetscFunctionBegin;
1641   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
1642   PetscFunctionReturn(0);
1643 }
1644 
1645 #undef __FUNCT__
1646 #define __FUNCT__ "DMLocalToGlobalEnd_Network"
1647 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
1648 {
1649   PetscErrorCode ierr;
1650   DM_Network     *network = (DM_Network*) dm->data;
1651 
1652   PetscFunctionBegin;
1653   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
1654   PetscFunctionReturn(0);
1655 }
1656