xref: /petsc/src/dm/impls/network/network.c (revision eab1376d19bf8559d4061dc648d5989f264e624e)
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 = PetscMalloc2(nleaves_sub,&ilocal_sub,nleaves_sub,&iremote_sub);CHKERRQ(ierr);
835   nleaves_sub = 0;
836   for (i = 0; i < nleaves; i++) {
837     if (ilocal_map[i] != -1) {
838       ilocal_sub[nleaves_sub] = ilocal_map[i];
839       iremote_sub[nleaves_sub] = iremote[i];
840       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
841       nleaves_sub += 1;
842     }
843   }
844   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
845   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
846 
847   /* Create new subSF */
848   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
849   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
850   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_OWN_POINTER);CHKERRQ(ierr);
851   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
852 
853   PetscFunctionReturn(0);
854 }
855 
856 
857 #undef __FUNCT__
858 #define __FUNCT__ "DMNetworkGetSupportingEdges"
859 /*@C
860   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
861 
862   Not Collective
863 
864   Input Parameters:
865 + dm - The DMNetwork object
866 - p  - the vertex point
867 
868   Output Paramters:
869 + nedges - number of edges connected to this vertex point
870 - edges  - List of edge points
871 
872   Level: intermediate
873 
874   Fortran Notes:
875   Since it returns an array, this routine is only available in Fortran 90, and you must
876   include petsc.h90 in your code.
877 
878 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
879 @*/
880 PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
881 {
882   PetscErrorCode ierr;
883   DM_Network     *network = (DM_Network*)dm->data;
884 
885   PetscFunctionBegin;
886   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
887   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
888   PetscFunctionReturn(0);
889 }
890 
891 #undef __FUNCT__
892 #define __FUNCT__ "DMNetworkGetConnectedNodes"
893 /*@C
894   DMNetworkGetConnectedNodes - Return the connected vertices for this edge point
895 
896   Not Collective
897 
898   Input Parameters:
899 + dm - The DMNetwork object
900 - p  - the edge point
901 
902   Output Paramters:
903 . vertices  - vertices connected to this edge
904 
905   Level: intermediate
906 
907   Fortran Notes:
908   Since it returns an array, this routine is only available in Fortran 90, and you must
909   include petsc.h90 in your code.
910 
911 .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
912 @*/
913 PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
914 {
915   PetscErrorCode ierr;
916   DM_Network     *network = (DM_Network*)dm->data;
917 
918   PetscFunctionBegin;
919   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
920   PetscFunctionReturn(0);
921 }
922 
923 #undef __FUNCT__
924 #define __FUNCT__ "DMNetworkIsGhostVertex"
925 /*@
926   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
927 
928   Not Collective
929 
930   Input Parameters:
931 + dm - The DMNetwork object
932 . p  - the vertex point
933 
934   Output Parameter:
935 . isghost - TRUE if the vertex is a ghost point
936 
937   Level: intermediate
938 
939 .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
940 @*/
941 PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
942 {
943   PetscErrorCode ierr;
944   DM_Network     *network = (DM_Network*)dm->data;
945   PetscInt       offsetg;
946   PetscSection   sectiong;
947 
948   PetscFunctionBegin;
949   *isghost = PETSC_FALSE;
950   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
951   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
952   if (offsetg < 0) *isghost = PETSC_TRUE;
953   PetscFunctionReturn(0);
954 }
955 
956 #undef __FUNCT__
957 #define __FUNCT__ "DMSetUp_Network"
958 PetscErrorCode DMSetUp_Network(DM dm)
959 {
960   PetscErrorCode ierr;
961   DM_Network     *network=(DM_Network*)dm->data;
962 
963   PetscFunctionBegin;
964   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
965   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
966 
967   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
968   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
969   PetscFunctionReturn(0);
970 }
971 
972 #undef __FUNCT__
973 #define __FUNCT__ "DMNetworkHasJacobian"
974 /*@
975     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
976                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
977 
978     Collective
979 
980     Input Parameters:
981 +   dm - The DMNetwork object
982 .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
983 -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
984 
985     Level: intermediate
986 
987 @*/
988 PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
989 {
990   DM_Network     *network=(DM_Network*)dm->data;
991 
992   PetscFunctionBegin;
993   network->userEdgeJacobian   = eflg;
994   network->userVertexJacobian = vflg;
995   PetscFunctionReturn(0);
996 }
997 
998 #undef __FUNCT__
999 #define __FUNCT__ "DMNetworkEdgeSetMatrix"
1000 /*@
1001     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
1002 
1003     Not Collective
1004 
1005     Input Parameters:
1006 +   dm - The DMNetwork object
1007 .   p  - the edge point
1008 -   J - array (size = 3) of Jacobian submatrices for this edge point:
1009         J[0]: this edge
1010         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes()
1011 
1012     Level: intermediate
1013 
1014 .seealso: DMNetworkVertexSetMatrix
1015 @*/
1016 PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
1017 {
1018   PetscErrorCode ierr;
1019   DM_Network     *network=(DM_Network*)dm->data;
1020 
1021   PetscFunctionBegin;
1022   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
1023   if (!network->Je) {
1024     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
1025   }
1026   network->Je[3*p]   = J[0];
1027   network->Je[3*p+1] = J[1];
1028   network->Je[3*p+2] = J[2];
1029   PetscFunctionReturn(0);
1030 }
1031 
1032 #undef __FUNCT__
1033 #define __FUNCT__ "DMNetworkVertexSetMatrix"
1034 /*@
1035     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
1036 
1037     Not Collective
1038 
1039     Input Parameters:
1040 +   dm - The DMNetwork object
1041 .   p  - the vertex point
1042 -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
1043         J[0]:       this vertex
1044         J[1+2*i]:   i-th supporting edge
1045         J[1+2*i+1]: i-th connected vertex
1046 
1047     Level: intermediate
1048 
1049 .seealso: DMNetworkEdgeSetMatrix
1050 @*/
1051 PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
1052 {
1053   PetscErrorCode ierr;
1054   DM_Network     *network=(DM_Network*)dm->data;
1055   PetscInt       i,*vptr,nedges,vStart,vEnd;
1056   const PetscInt *edges;
1057 
1058   PetscFunctionBegin;
1059   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1060 
1061   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1062 
1063   if (!network->Jv) {
1064     PetscInt nNodes = network->nNodes,nedges_total;
1065     /* count nvertex_total */
1066     nedges_total = 0;
1067     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1068     ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr);
1069 
1070     vptr[0] = 0;
1071     for (i=0; i<nNodes; i++) {
1072       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1073       nedges_total += nedges;
1074       vptr[i+1] = vptr[i] + 2*nedges + 1;
1075     }
1076 
1077     ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr);
1078     network->Jvptr = vptr;
1079   }
1080 
1081   vptr = network->Jvptr;
1082   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
1083 
1084   /* Set Jacobian for each supporting edge and connected vertex */
1085   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1086   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1087   PetscFunctionReturn(0);
1088 }
1089 
1090 #undef __FUNCT__
1091 #define __FUNCT__ "MatSetPreallocationDenseblock_private"
1092 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1093 {
1094   PetscErrorCode ierr;
1095   PetscInt       j;
1096   PetscScalar    val=(PetscScalar)ncols;
1097 
1098   PetscFunctionBegin;
1099   if (!ghost) {
1100     for (j=0; j<nrows; j++) {
1101       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1102     }
1103   } else {
1104     for (j=0; j<nrows; j++) {
1105       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1106     }
1107   }
1108   PetscFunctionReturn(0);
1109 }
1110 
1111 #undef __FUNCT__
1112 #define __FUNCT__ "MatSetPreallocationUserblock_private"
1113 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1114 {
1115   PetscErrorCode ierr;
1116   PetscInt       j,ncols_u;
1117   PetscScalar    val;
1118 
1119   PetscFunctionBegin;
1120   if (!ghost) {
1121     for (j=0; j<nrows; j++) {
1122       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1123       val = (PetscScalar)ncols_u;
1124       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1125       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1126     }
1127   } else {
1128     for (j=0; j<nrows; j++) {
1129       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1130       val = (PetscScalar)ncols_u;
1131       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1132       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1133     }
1134   }
1135   PetscFunctionReturn(0);
1136 }
1137 
1138 #undef __FUNCT__
1139 #define __FUNCT__ "MatSetPreallocationblock_private"
1140 PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
1141 {
1142   PetscErrorCode ierr;
1143 
1144   PetscFunctionBegin;
1145   if (Ju) {
1146     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
1147   } else {
1148     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
1149   }
1150   PetscFunctionReturn(0);
1151 }
1152 
1153 #undef __FUNCT__
1154 #define __FUNCT__ "MatSetDenseblock_private"
1155 PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1156 {
1157   PetscErrorCode ierr;
1158   PetscInt       j,*cols;
1159   PetscScalar    *zeros;
1160 
1161   PetscFunctionBegin;
1162   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1163   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1164   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1165   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
1166   PetscFunctionReturn(0);
1167 }
1168 
1169 #undef __FUNCT__
1170 #define __FUNCT__ "MatSetUserblock_private"
1171 PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1172 {
1173   PetscErrorCode ierr;
1174   PetscInt       j,M,N,row,col,ncols_u;
1175   const PetscInt *cols;
1176   PetscScalar    zero=0.0;
1177 
1178   PetscFunctionBegin;
1179   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
1180   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);
1181 
1182   for (row=0; row<nrows; row++) {
1183     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
1184     for (j=0; j<ncols_u; j++) {
1185       col = cols[j] + cstart;
1186       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
1187     }
1188     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
1189   }
1190   PetscFunctionReturn(0);
1191 }
1192 
1193 #undef __FUNCT__
1194 #define __FUNCT__ "MatSetblock_private"
1195 PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1196 {
1197   PetscErrorCode ierr;
1198 
1199   PetscFunctionBegin;
1200   if (Ju) {
1201     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1202   } else {
1203     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1204   }
1205   PetscFunctionReturn(0);
1206 }
1207 
1208 
1209 #undef __FUNCT__
1210 #define __FUNCT__ "CreateSubGlobalToLocalMapping_private"
1211 /* 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.
1212 */
1213 PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
1214 {
1215   PetscErrorCode ierr;
1216   PetscInt       i, size, dof;
1217   PetscInt       *glob2loc;
1218 
1219   PetscFunctionBegin;
1220   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
1221   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
1222 
1223   for (i = 0; i < size; i++) {
1224     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
1225     dof = (dof >= 0) ? dof : -(dof + 1);
1226     glob2loc[i] = dof;
1227   }
1228 
1229   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
1230 #if 0
1231   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
1232 #endif
1233   PetscFunctionReturn(0);
1234 }
1235 
1236 #undef __FUNCT__
1237 #define __FUNCT__ "DMCreateMatrix_Network"
1238 PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
1239 {
1240   PetscErrorCode ierr;
1241   PetscMPIInt    rank, size;
1242   DM_Network     *network = (DM_Network*) dm->data;
1243   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1244   PetscInt       cstart,ncols,j,e,v;
1245   PetscBool      ghost,ghost_vc,ghost2,isNest;
1246   Mat            Juser;
1247   PetscSection   sectionGlobal;
1248   PetscInt       nedges,*vptr=NULL,vc,*rows_v;
1249   const PetscInt *edges,*cone;
1250   MPI_Comm       comm;
1251   MatType        mtype;
1252   Vec            vd_nz,vo_nz;
1253   PetscInt       *dnnz,*onnz;
1254   PetscScalar    *vdnz,*vonz;
1255 
1256   PetscFunctionBegin;
1257   mtype = dm->mattype;
1258   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
1259 
1260   if (isNest) {
1261     /* ierr = DMCreateMatrix_Network_Nest(); */
1262     PetscInt   eDof, vDof;
1263     Mat        j11, j12, j21, j22, bA[2][2];
1264     ISLocalToGlobalMapping eISMap, vISMap;
1265 
1266     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1267     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
1268     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
1269 
1270     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
1271     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
1272 
1273 #if 0
1274     printf("rank[%d] vdof = %d. edof = %d\n", rank, vDof, eDof);
1275 #endif
1276 
1277     ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr);
1278     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1279     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
1280 
1281     ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr);
1282     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
1283     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
1284 
1285     ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr);
1286     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1287     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
1288 
1289     ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr);
1290     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
1291     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
1292 
1293     bA[0][0] = j11;
1294     bA[0][1] = j12;
1295     bA[1][0] = j21;
1296     bA[1][1] = j22;
1297 
1298     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
1299     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
1300 
1301     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
1302     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
1303     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
1304     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
1305 
1306     ierr = MatSetUp(j11);CHKERRQ(ierr);
1307     ierr = MatSetUp(j12);CHKERRQ(ierr);
1308     ierr = MatSetUp(j21);CHKERRQ(ierr);
1309     ierr = MatSetUp(j22);CHKERRQ(ierr);
1310 
1311     ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
1312     ierr = MatSetUp(*J);CHKERRQ(ierr);
1313     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
1314     ierr = MatDestroy(&j11);CHKERRQ(ierr);
1315     ierr = MatDestroy(&j12);CHKERRQ(ierr);
1316     ierr = MatDestroy(&j21);CHKERRQ(ierr);
1317     ierr = MatDestroy(&j22);CHKERRQ(ierr);
1318 
1319     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1320     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1321 
1322     /* Free structures */
1323     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
1324     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
1325 
1326     PetscFunctionReturn(0);
1327   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1328     /* user does not provide Jacobian blocks */
1329     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1330     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1331     PetscFunctionReturn(0);
1332   }
1333 
1334   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
1335   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1336   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1337   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
1338 
1339   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
1340   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
1341 
1342   /* (1) Set matrix preallocation */
1343   /*------------------------------*/
1344   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1345   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1346   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1347   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1348   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1349   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1350 
1351   /* Set preallocation for edges */
1352   /*-----------------------------*/
1353   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1354 
1355   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1356   for (e=eStart; e<eEnd; e++) {
1357     /* Get row indices */
1358     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1359     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1360     if (nrows) {
1361       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1362       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1363 
1364       /* Set preallocation for conntected vertices */
1365       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1366       for (v=0; v<2; v++) {
1367         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1368 
1369         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1370         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
1371         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1372       }
1373 
1374       /* Set preallocation for edge self */
1375       cstart = rstart;
1376       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1377       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1378     }
1379   }
1380 
1381   /* Set preallocation for vertices */
1382   /*--------------------------------*/
1383   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1384   if (vEnd - vStart) {
1385     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1386     vptr = network->Jvptr;
1387     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1388   }
1389 
1390   for (v=vStart; v<vEnd; v++) {
1391     /* Get row indices */
1392     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1393     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1394     if (!nrows) continue;
1395 
1396     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1397     if (ghost) {
1398       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1399     } else {
1400       rows_v = rows;
1401     }
1402 
1403     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1404 
1405     /* Get supporting edges and connected vertices */
1406     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1407 
1408     for (e=0; e<nedges; e++) {
1409       /* Supporting edges */
1410       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1411       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1412 
1413       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1414       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1415 
1416       /* Connected vertices */
1417       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1418       vc = (v == cone[0]) ? cone[1]:cone[0];
1419       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1420 
1421       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1422 
1423       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1424       if (ghost_vc||ghost) {
1425         ghost2 = PETSC_TRUE;
1426       } else {
1427         ghost2 = PETSC_FALSE;
1428       }
1429       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1430     }
1431 
1432     /* Set preallocation for vertex self */
1433     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1434     if (!ghost) {
1435       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1436       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1437       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1438     }
1439     if (ghost) {
1440       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1441     }
1442   }
1443 
1444   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1445   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
1446 
1447   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
1448 
1449   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1450   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1451 
1452   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1453   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1454   for (j=0; j<localSize; j++) {
1455     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1456     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1457   }
1458   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1459   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1460   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1461   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1462 
1463   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
1464   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
1465   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
1466 
1467   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
1468 
1469   /* (2) Set matrix entries for edges */
1470   /*----------------------------------*/
1471   for (e=eStart; e<eEnd; e++) {
1472     /* Get row indices */
1473     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1474     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1475     if (nrows) {
1476       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1477 
1478       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1479 
1480       /* Set matrix entries for conntected vertices */
1481       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1482       for (v=0; v<2; v++) {
1483         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1484         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1485 
1486         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1487         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1488       }
1489 
1490       /* Set matrix entries for edge self */
1491       cstart = rstart;
1492       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1493       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1494     }
1495   }
1496 
1497   /* Set matrix entries for vertices */
1498   /*---------------------------------*/
1499   for (v=vStart; v<vEnd; v++) {
1500     /* Get row indices */
1501     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1502     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1503     if (!nrows) continue;
1504 
1505     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1506     if (ghost) {
1507       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1508     } else {
1509       rows_v = rows;
1510     }
1511     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1512 
1513     /* Get supporting edges and connected vertices */
1514     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1515 
1516     for (e=0; e<nedges; e++) {
1517       /* Supporting edges */
1518       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1519       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1520 
1521       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1522       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1523 
1524       /* Connected vertices */
1525       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1526       vc = (v == cone[0]) ? cone[1]:cone[0];
1527 
1528       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
1529       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1530 
1531       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1532       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1533     }
1534 
1535     /* Set matrix entries for vertex self */
1536     if (!ghost) {
1537       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1538       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1539       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1540     }
1541     if (ghost) {
1542       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1543     }
1544   }
1545   ierr = PetscFree(rows);CHKERRQ(ierr);
1546 
1547   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1548   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1549 
1550   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
1551   PetscFunctionReturn(0);
1552 }
1553 
1554 #undef __FUNCT__
1555 #define __FUNCT__ "DMDestroy_Network"
1556 PetscErrorCode DMDestroy_Network(DM dm)
1557 {
1558   PetscErrorCode ierr;
1559   DM_Network     *network = (DM_Network*) dm->data;
1560 
1561   PetscFunctionBegin;
1562   if (--network->refct > 0) PetscFunctionReturn(0);
1563   if (network->Je) {
1564     ierr = PetscFree(network->Je);CHKERRQ(ierr);
1565   }
1566   if (network->Jv) {
1567     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
1568     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
1569   }
1570   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
1571   network->edges = NULL;
1572   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
1573   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
1574 
1575   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
1576   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
1577   ierr = PetscFree(network->header);CHKERRQ(ierr);
1578   ierr = PetscFree(network);CHKERRQ(ierr);
1579   PetscFunctionReturn(0);
1580 }
1581 
1582 #undef __FUNCT__
1583 #define __FUNCT__ "DMView_Network"
1584 PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
1585 {
1586   PetscErrorCode ierr;
1587   DM_Network     *network = (DM_Network*) dm->data;
1588 
1589   PetscFunctionBegin;
1590   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
1591   PetscFunctionReturn(0);
1592 }
1593 
1594 #undef __FUNCT__
1595 #define __FUNCT__ "DMGlobalToLocalBegin_Network"
1596 PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
1597 {
1598   PetscErrorCode ierr;
1599   DM_Network     *network = (DM_Network*) dm->data;
1600 
1601   PetscFunctionBegin;
1602   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
1603   PetscFunctionReturn(0);
1604 }
1605 
1606 #undef __FUNCT__
1607 #define __FUNCT__ "DMGlobalToLocalEnd_Network"
1608 PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
1609 {
1610   PetscErrorCode ierr;
1611   DM_Network     *network = (DM_Network*) dm->data;
1612 
1613   PetscFunctionBegin;
1614   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
1615   PetscFunctionReturn(0);
1616 }
1617 
1618 #undef __FUNCT__
1619 #define __FUNCT__ "DMLocalToGlobalBegin_Network"
1620 PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
1621 {
1622   PetscErrorCode ierr;
1623   DM_Network     *network = (DM_Network*) dm->data;
1624 
1625   PetscFunctionBegin;
1626   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
1627   PetscFunctionReturn(0);
1628 }
1629 
1630 #undef __FUNCT__
1631 #define __FUNCT__ "DMLocalToGlobalEnd_Network"
1632 PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
1633 {
1634   PetscErrorCode ierr;
1635   DM_Network     *network = (DM_Network*) dm->data;
1636 
1637   PetscFunctionBegin;
1638   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
1639   PetscFunctionReturn(0);
1640 }
1641