Actual source code: network.c

petsc-3.5.4 2015-05-23
Report Typos and Errors
  1: #include <petsc-private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
  2: #include <petscdmplex.h>
  3: #include <petscsf.h>

  7: /*@
  8:   DMNetworkSetSizes - Sets the local and global vertices and edges.

 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)

 19:    Notes
 20:    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.

 22:    You cannot change the sizes once they have been set

 24:    Level: intermediate

 26: .seealso: DMNetworkCreate
 27: @*/
 28: PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE)
 29: {
 31:   DM_Network     *network = (DM_Network*) dm->data;
 32:   PetscInt       a[2],b[2];

 38:   if (NV > 0 && nV > NV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local vertex size %D cannot be larger than global vertex size %D",nV,NV);
 39:   if (NE > 0 && nE > NE) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local edge size %D cannot be larger than global edge size %D",nE,NE);
 40:   if ((network->nNodes >= 0 || network->NNodes >= 0) && (network->nNodes != nV || network->NNodes != NV)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset vertex sizes to %D local %D global after previously setting them to %D local %D global",nV,NV,network->nNodes,network->NNodes);
 41:   if ((network->nEdges >= 0 || network->NEdges >= 0) && (network->nEdges != nE || network->NEdges != NE)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset edge sizes to %D local %D global after previously setting them to %D local %D global",nE,NE,network->nEdges,network->NEdges);
 42:   if (NE < 0 || NV < 0) {
 43:     a[0] = nV; a[1] = nE;
 44:     MPI_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));
 45:     NV = b[0]; NE = b[1];
 46:   }
 47:   network->nNodes = nV;
 48:   network->NNodes = NV;
 49:   network->nEdges = nE;
 50:   network->NEdges = NE;
 51:   return(0);
 52: }

 56: /*@
 57:   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network

 59:   Logically collective on DM

 61:   Input Parameters:
 62: . edges - list of edges

 64:   Notes:
 65:   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
 66:   not be destroyed before the call to DMNetworkLayoutSetUp

 68:   Level: intermediate

 70: .seealso: DMNetworkCreate, DMNetworkSetSizes
 71: @*/
 72: PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[])
 73: {
 74:   DM_Network *network = (DM_Network*) dm->data;
 75: 
 77:   network->edges = edgelist;
 78:   return(0);
 79: }

 83: /*@
 84:   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network

 86:   Collective on DM

 88:   Input Parameters
 89: . DM - the dmnetwork object

 91:   Notes:
 92:   This routine should be called after the network sizes and edgelists have been provided. It creates
 93:   the bare layout of the network and sets up the network to begin insertion of components.

 95:   All the components should be registered before calling this routine.

 97:   Level: intermediate

 99: .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
100: @*/
101: PetscErrorCode DMNetworkLayoutSetUp(DM dm)
102: {
104:   DM_Network     *network = (DM_Network*) dm->data;
105:   PetscInt       dim = 1; /* One dimensional network */
106:   PetscInt       numCorners=2;
107:   PetscInt       spacedim=2;
108:   double         *vertexcoords=NULL;
109:   PetscInt       i;
110:   PetscInt       ndata;

113:   if (network->nNodes) {
114:     PetscMalloc1(numCorners*network->nNodes,&vertexcoords);
115:   }
116:   DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);
117:   if (network->nNodes) {
118:     PetscFree(vertexcoords);
119:   }
120:   DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);
121:   DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);
122:   DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);
123: 
124:   PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);
125:   PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);
126:   PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);
127:   PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);

129:   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
130:   PetscMalloc1((network->pEnd-network->pStart),&network->header);
131:   for (i = network->pStart; i < network->pEnd; i++) {
132:     network->header[i].ndata = 0;
133:     ndata = network->header[i].ndata;
134:     PetscSectionAddDof(network->DataSection,i,network->dataheadersize);
135:     network->header[i].offset[ndata] = 0;
136:   }
137:   PetscMalloc1((network->pEnd-network->pStart),&network->cvalue);
138:   return(0);
139: }

143: /*@
144:   DMNetworkRegisterComponent - Registers the network component

146:   Logically collective on DM

148:   Input Parameters
149: + dm   - the network object
150: . name - the component name
151: - size - the storage size in bytes for this component data

153:    Output Parameters
154: .   key - an integer key that defines the component

156:    Notes
157:    This routine should be called by all processors before calling DMNetworkLayoutSetup().

159:    Level: intermediate

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;


173:   for (i=0; i < network->ncomponent; i++) {
174:     PetscStrcmp(component->name,name,&flg);
175:     if (flg) {
176:       *key = i;
177:       return(0);
178:     }
179:   }
180: 
181:   PetscStrcpy(component->name,name);
182:   component->size = size/sizeof(DMNetworkComponentGenericDataType);
183:   *key = network->ncomponent;
184:   network->ncomponent++;
185:   return(0);
186: }

190: /*@
191:   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.

193:   Not Collective

195:   Input Parameters:
196: + dm - The DMNetwork object

198:   Output Paramters:
199: + vStart - The first vertex point
200: - vEnd   - One beyond the last vertex point

202:   Level: intermediate

204: .seealso: DMNetworkGetEdgeRange
205: @*/
206: PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
207: {
208:   DM_Network     *network = (DM_Network*)dm->data;

211:   if (vStart) *vStart = network->vStart;
212:   if (vEnd) *vEnd = network->vEnd;
213:   return(0);
214: }

218: /*@
219:   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.

221:   Not Collective

223:   Input Parameters:
224: + dm - The DMNetwork object

226:   Output Paramters:
227: + eStart - The first edge point
228: - eEnd   - One beyond the last edge point

230:   Level: intermediate

232: .seealso: DMNetworkGetVertexRange
233: @*/
234: PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
235: {
236:   DM_Network     *network = (DM_Network*)dm->data;

239:   if (eStart) *eStart = network->eStart;
240:   if (eEnd) *eEnd = network->eEnd;
241:   return(0);
242: }

246: /*@
247:   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)

249:   Not Collective

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

257:   Level: intermediate

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: 
270:   header->size[header->ndata] = component.size;
271:   PetscSectionAddDof(network->DataSection,p,component.size);
272:   header->key[header->ndata] = componentkey;
273:   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];

275:   cvalue->data[header->ndata] = (void*)compvalue;
276:   header->ndata++;
277:   return(0);
278: }

282: /*@
283:   DMNetworkGetNumComponents - Get the number of components at a vertex/edge

285:   Not Collective 

287:   Input Parameters:
288: + dm - The DMNetwork object
289: . p  - vertex/edge point

291:   Output Parameters:
292: . numcomponents - Number of components at the vertex/edge

294:   Level: intermediate

296: .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
297: @*/
298: PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
299: {
301:   PetscInt       offset;
302:   DM_Network     *network = (DM_Network*)dm->data;

305:   PetscSectionGetOffset(network->DataSection,p,&offset);
306:   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
307:   return(0);
308: }

312: /*@
313:   DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the 
314:                                     component value from the component data array

316:   Not Collective

318:   Input Parameters:
319: + dm      - The DMNetwork object
320: . p       - vertex/edge point
321: - compnum - component number
322:         
323:   Output Parameters:
324: + compkey - the key obtained when registering the component
325: - offset  - offset into the component data array associated with the vertex/edge point

327:   Notes:
328:   Typical usage:

330:   DMNetworkGetComponentDataArray(dm, &arr);
331:   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
332:   Loop over vertices or edges
333:     DMNetworkGetNumComponents(dm,v,&numcomps);
334:     Loop over numcomps
335:       DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset);
336:       compdata = (UserCompDataType)(arr+offset);
337:   
338:   Level: intermediate

340: .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 
341: @*/
342: PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
343: {
345:   PetscInt       offsetp;
346:   DMNetworkComponentHeader header;
347:   DM_Network     *network = (DM_Network*)dm->data;

350:   PetscSectionGetOffset(network->DataSection,p,&offsetp);
351:   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
352:   *compkey = header->key[compnum];
353:   *offset  = offsetp+network->dataheadersize+header->offset[compnum];
354:   return(0);
355: }

359: /*@
360:   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.

362:   Not Collective

364:   Input Parameters:
365: + dm     - The DMNetwork object
366: - p      - the edge/vertex point

368:   Output Parameters:
369: . offset - the offset

371:   Level: intermediate

373: .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
374: @*/
375: PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
376: {
378:   DM_Network     *network = (DM_Network*)dm->data;

381:   PetscSectionGetOffset(network->DofSection,p,offset);
382:   return(0);
383: }

387: /*@
388:   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.

390:   Not Collective

392:   Input Parameters:
393: + dm      - The DMNetwork object
394: - p       - the edge/vertex point

396:   Output Parameters:
397: . offsetg - the offset

399:   Level: intermediate

401: .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
402: @*/
403: PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
404: {
406:   DM_Network     *network = (DM_Network*)dm->data;

409:   PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);
410:   return(0);
411: }

415: /*@ 
416:   DMNetworkAddNumVariables - Add number of variables associated with a given point.

418:   Not Collective

420:   Input Parameters:
421: + dm   - The DMNetworkObject
422: . p    - the vertex/edge point
423: - nvar - number of additional variables

425:   Level: intermediate

427: .seealso: DMNetworkSetNumVariables
428: @*/
429: PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
430: {
432:   DM_Network     *network = (DM_Network*)dm->data;

435:   PetscSectionAddDof(network->DofSection,p,nvar);
436:   return(0);
437: }

441: /*@ 
442:   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.

444:   Not Collective

446:   Input Parameters:
447: + dm   - The DMNetworkObject
448: . p    - the vertex/edge point
449: - nvar - number of variables

451:   Level: intermediate

453: .seealso: DMNetworkAddNumVariables
454: @*/
455: PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
456: {
458:   DM_Network     *network = (DM_Network*)dm->data;

461:   PetscSectionSetDof(network->DofSection,p,nvar);
462:   return(0);
463: }

465: /* Sets up the array that holds the data for all components and its associated section. This
466:    function is called during DMSetUp() */
469: PetscErrorCode DMNetworkComponentSetUp(DM dm)
470: {
471:   PetscErrorCode              ierr;
472:   DM_Network     *network = (DM_Network*)dm->data;
473:   PetscInt                    arr_size;
474:   PetscInt                    p,offset,offsetp;
475:   DMNetworkComponentHeader header;
476:   DMNetworkComponentValue  cvalue;
477:   DMNetworkComponentGenericDataType      *componentdataarray;
478:   PetscInt ncomp, i;

481:   PetscSectionSetUp(network->DataSection);
482:   PetscSectionGetStorageSize(network->DataSection,&arr_size);
483:   PetscMalloc1(arr_size,&network->componentdataarray);
484:   componentdataarray = network->componentdataarray;
485:   for (p = network->pStart; p < network->pEnd; p++) {
486:     PetscSectionGetOffset(network->DataSection,p,&offsetp);
487:     /* Copy header */
488:     header = &network->header[p];
489:     PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));
490:     /* Copy data */
491:     cvalue = &network->cvalue[p];
492:     ncomp = header->ndata;
493:     for (i = 0; i < ncomp; i++) {
494:       offset = offsetp + network->dataheadersize + header->offset[i];
495:       PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));
496:     }
497:   }
498:   return(0);
499: }

501: /* Sets up the section for dofs. This routine is called during DMSetUp() */
504: PetscErrorCode DMNetworkVariablesSetUp(DM dm)
505: {
507:   DM_Network     *network = (DM_Network*)dm->data;

510:   PetscSectionSetUp(network->DofSection);
511:   return(0);
512: }

516: /*@C
517:   DMNetworkGetComponentDataArray - Returns the component data array

519:   Not Collective

521:   Input Parameters:
522: . dm - The DMNetwork Object

524:   Output Parameters:
525: . componentdataarray - array that holds data for all components

527:   Level: intermediate

529: .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents
530: @*/
531: PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
532: {
533:   DM_Network     *network = (DM_Network*)dm->data;

536:   *componentdataarray = network->componentdataarray;
537:   return(0);
538: }

542: /*@
543:   DMNetworkDistribute - Distributes the network and moves associated component data.

545:   Collective

547:   Input Parameter:
548: + oldDM - the original DMNetwork object
549: . partitioner - The partitioning package, or NULL for the default
550: - overlap - The overlap of partitions, 0 is the default

552:   Output Parameter:
553: . distDM - the distributed DMNetwork object

555:   Notes:
556:   This routine should be called only when using multiple processors.

558:   Distributes the network with a non-overlapping partitioning of the edges.

560:   Level: intermediate

562: .seealso: DMNetworkCreate
563: @*/
564: PetscErrorCode DMNetworkDistribute(DM oldDM, const char partitioner[], PetscInt overlap,DM *distDM)
565: {
567:   DM_Network     *oldDMnetwork = (DM_Network*)oldDM->data;
568:   PetscSF        pointsf;
569:   DM             newDM;
570:   DM_Network     *newDMnetwork;

573:   DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);
574:   newDMnetwork = (DM_Network*)newDM->data;
575:   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
576:   /* Distribute plex dm and dof section */
577:   DMPlexDistribute(oldDMnetwork->plex,partitioner,overlap,&pointsf,&newDMnetwork->plex);
578:   /* Distribute dof section */
579:   PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);
580:   PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);
581:   PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);
582:   /* Distribute data and associated section */
583:   DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPI_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);
584:   /* Destroy point SF */
585:   PetscSFDestroy(&pointsf);
586: 
587:   PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);
588:   DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);
589:   DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);
590:   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
591:   newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart;
592:   newDMnetwork->NNodes = oldDMnetwork->NNodes;
593:   newDMnetwork->NEdges = oldDMnetwork->NEdges;
594:   /* Set Dof section as the default section for dm */
595:   DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);
596:   DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);

598:   *distDM = newDM;
599:   return(0);
600: }

604: /*@C
605:   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point

607:   Not Collective

609:   Input Parameters:
610: + dm - The DMNetwork object
611: - p  - the vertex point

613:   Output Paramters:
614: + nedges - number of edges connected to this vertex point
615: - edges  - List of edge points

617:   Level: intermediate

619:   Fortran Notes:
620:   Since it returns an array, this routine is only available in Fortran 90, and you must
621:   include petsc.h90 in your code.

623: .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
624: @*/
625: PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
626: {
628:   DM_Network     *network = (DM_Network*)dm->data;

631:   DMPlexGetSupportSize(network->plex,vertex,nedges);
632:   DMPlexGetSupport(network->plex,vertex,edges);
633:   return(0);
634: }

638: /*@C
639:   DMNetworkGetConnectedNodes - Return the connected edges for this edge point

641:   Not Collective

643:   Input Parameters:
644: + dm - The DMNetwork object
645: - p  - the edge point

647:   Output Paramters:
648: . vertices  - vertices connected to this edge

650:   Level: intermediate

652:   Fortran Notes:
653:   Since it returns an array, this routine is only available in Fortran 90, and you must
654:   include petsc.h90 in your code.

656: .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
657: @*/
658: PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
659: {
661:   DM_Network     *network = (DM_Network*)dm->data;

664:   DMPlexGetCone(network->plex,edge,vertices);
665:   return(0);
666: }

670: /*@
671:   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex

673:   Not Collective

675:   Input Parameters:
676: + dm - The DMNetwork object
677: . p  - the vertex point

679:   Output Parameter:
680: . isghost - TRUE if the vertex is a ghost point 

682:   Level: intermediate

684: .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
685: @*/
686: PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
687: {
689:   DM_Network     *network = (DM_Network*)dm->data;
690:   PetscInt       offsetg;
691:   PetscSection   sectiong;

694:   *isghost = PETSC_FALSE;
695:   DMGetDefaultGlobalSection(network->plex,&sectiong);
696:   PetscSectionGetOffset(sectiong,p,&offsetg);
697:   if (offsetg < 0) *isghost = PETSC_TRUE;
698:   return(0);
699: }

703: PetscErrorCode DMSetUp_Network(DM dm)
704: {
706:   DM_Network     *network=(DM_Network*)dm->data;

709:   DMNetworkComponentSetUp(dm);
710:   DMNetworkVariablesSetUp(dm);

712:   DMSetDefaultSection(network->plex,network->DofSection);
713:   DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);
714:   return(0);
715: }

719: PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
720: {
722:   DM_Network     *network = (DM_Network*) dm->data;

725:   DMCreateMatrix(network->plex,J);
726:   MatSetDM(*J,dm);
727:   return(0);
728: }

732: PetscErrorCode DMDestroy_Network(DM dm)
733: {
735:   DM_Network     *network = (DM_Network*) dm->data;

738:   DMDestroy(&network->plex);
739:   network->edges = NULL;
740:   PetscSectionDestroy(&network->DataSection);
741:   PetscSectionDestroy(&network->DofSection);
742:   /*  PetscSectionDestroy(&network->GlobalDofSection); */
743:   PetscFree(network->componentdataarray);
744:   PetscFree(network->cvalue);
745:   PetscFree(network->header);
746:   PetscFree(network);
747:   return(0);
748: }

752: PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
753: {
755:   DM_Network     *network = (DM_Network*) dm->data;

758:   DMView(network->plex,viewer);
759:   return(0);
760: }

764: PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
765: {
767:   DM_Network     *network = (DM_Network*) dm->data;

770:   DMGlobalToLocalBegin(network->plex,g,mode,l);
771:   return(0);
772: }

776: PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
777: {
779:   DM_Network     *network = (DM_Network*) dm->data;

782:   DMGlobalToLocalEnd(network->plex,g,mode,l);
783:   return(0);
784: }

788: PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
789: {
791:   DM_Network     *network = (DM_Network*) dm->data;

794:   DMLocalToGlobalBegin(network->plex,l,mode,g);
795:   return(0);
796: }

800: PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
801: {
803:   DM_Network     *network = (DM_Network*) dm->data;

806:   DMLocalToGlobalEnd(network->plex,l,mode,g);
807:   return(0);
808: }