xref: /petsc/src/dm/impls/network/network.c (revision 83b2e829dcb382aebd416aff29be07233908db19)
1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h>  /*I  "petscdmnetwork.h"  I*/
25f2c45f1SShri Abhyankar #include <petscdmplex.h>
35f2c45f1SShri Abhyankar #include <petscsf.h>
45f2c45f1SShri Abhyankar 
55f2c45f1SShri Abhyankar #undef __FUNCT__
65f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkSetSizes"
75f2c45f1SShri Abhyankar /*@
85f2c45f1SShri Abhyankar   DMNetworkSetSizes - Sets the local and global vertices and edges.
95f2c45f1SShri Abhyankar 
105f2c45f1SShri Abhyankar   Collective on DM
115f2c45f1SShri Abhyankar 
125f2c45f1SShri Abhyankar   Input Parameters:
135f2c45f1SShri Abhyankar + dm - the dm object
145f2c45f1SShri Abhyankar . nV - number of local vertices
155f2c45f1SShri Abhyankar . nE - number of local edges
165f2c45f1SShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE)
175f2c45f1SShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE)
185f2c45f1SShri Abhyankar 
195f2c45f1SShri Abhyankar    Notes
205f2c45f1SShri Abhyankar    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
215f2c45f1SShri Abhyankar 
225f2c45f1SShri Abhyankar    You cannot change the sizes once they have been set
235f2c45f1SShri Abhyankar 
245f2c45f1SShri Abhyankar    Level: intermediate
255f2c45f1SShri Abhyankar 
265f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
275f2c45f1SShri Abhyankar @*/
285f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE)
295f2c45f1SShri Abhyankar {
305f2c45f1SShri Abhyankar   PetscErrorCode ierr;
315f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
325f2c45f1SShri Abhyankar   PetscInt       a[2],b[2];
335f2c45f1SShri Abhyankar 
345f2c45f1SShri Abhyankar   PetscFunctionBegin;
355f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
365f2c45f1SShri Abhyankar   if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4);
375f2c45f1SShri Abhyankar   if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5);
385f2c45f1SShri Abhyankar   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);
395f2c45f1SShri Abhyankar   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);
405f2c45f1SShri Abhyankar   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);
415f2c45f1SShri Abhyankar   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);
425f2c45f1SShri Abhyankar   if (NE < 0 || NV < 0) {
435f2c45f1SShri Abhyankar     a[0] = nV; a[1] = nE;
44b2566f29SBarry Smith     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
455f2c45f1SShri Abhyankar     NV = b[0]; NE = b[1];
465f2c45f1SShri Abhyankar   }
475f2c45f1SShri Abhyankar   network->nNodes = nV;
485f2c45f1SShri Abhyankar   network->NNodes = NV;
495f2c45f1SShri Abhyankar   network->nEdges = nE;
505f2c45f1SShri Abhyankar   network->NEdges = NE;
515f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
525f2c45f1SShri Abhyankar }
535f2c45f1SShri Abhyankar 
545f2c45f1SShri Abhyankar #undef __FUNCT__
555f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkSetEdgeList"
565f2c45f1SShri Abhyankar /*@
575f2c45f1SShri Abhyankar   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
585f2c45f1SShri Abhyankar 
595f2c45f1SShri Abhyankar   Logically collective on DM
605f2c45f1SShri Abhyankar 
615f2c45f1SShri Abhyankar   Input Parameters:
625f2c45f1SShri Abhyankar . edges - list of edges
635f2c45f1SShri Abhyankar 
645f2c45f1SShri Abhyankar   Notes:
655f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
665f2c45f1SShri Abhyankar   not be destroyed before the call to DMNetworkLayoutSetUp
675f2c45f1SShri Abhyankar 
685f2c45f1SShri Abhyankar   Level: intermediate
695f2c45f1SShri Abhyankar 
705f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes
715f2c45f1SShri Abhyankar @*/
725f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[])
735f2c45f1SShri Abhyankar {
745f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
755f2c45f1SShri Abhyankar 
765f2c45f1SShri Abhyankar   PetscFunctionBegin;
775f2c45f1SShri Abhyankar   network->edges = edgelist;
785f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
795f2c45f1SShri Abhyankar }
805f2c45f1SShri Abhyankar 
815f2c45f1SShri Abhyankar #undef __FUNCT__
825f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkLayoutSetUp"
835f2c45f1SShri Abhyankar /*@
845f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
855f2c45f1SShri Abhyankar 
865f2c45f1SShri Abhyankar   Collective on DM
875f2c45f1SShri Abhyankar 
885f2c45f1SShri Abhyankar   Input Parameters
895f2c45f1SShri Abhyankar . DM - the dmnetwork object
905f2c45f1SShri Abhyankar 
915f2c45f1SShri Abhyankar   Notes:
925f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
935f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
945f2c45f1SShri Abhyankar 
955f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
965f2c45f1SShri Abhyankar 
975f2c45f1SShri Abhyankar   Level: intermediate
985f2c45f1SShri Abhyankar 
995f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
1005f2c45f1SShri Abhyankar @*/
1015f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm)
1025f2c45f1SShri Abhyankar {
1035f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1045f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
1055f2c45f1SShri Abhyankar   PetscInt       dim = 1; /* One dimensional network */
1065f2c45f1SShri Abhyankar   PetscInt       numCorners=2;
1075f2c45f1SShri Abhyankar   PetscInt       spacedim=2;
1085f2c45f1SShri Abhyankar   double         *vertexcoords=NULL;
1095f2c45f1SShri Abhyankar   PetscInt       i;
1105f2c45f1SShri Abhyankar   PetscInt       ndata;
1115f2c45f1SShri Abhyankar 
1125f2c45f1SShri Abhyankar   PetscFunctionBegin;
1135f2c45f1SShri Abhyankar   if (network->nNodes) {
1149045477aSSatish Balay     ierr = PetscCalloc1(numCorners*network->nNodes,&vertexcoords);CHKERRQ(ierr);
1155f2c45f1SShri Abhyankar   }
1165f2c45f1SShri Abhyankar   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
1175f2c45f1SShri Abhyankar   if (network->nNodes) {
1185f2c45f1SShri Abhyankar     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
1195f2c45f1SShri Abhyankar   }
1205f2c45f1SShri Abhyankar   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
1215f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
1225f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
1235f2c45f1SShri Abhyankar 
1245f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
1255f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
1265f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
1275f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
1285f2c45f1SShri Abhyankar 
1295f2c45f1SShri Abhyankar   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
1306caa05f4SBarry Smith   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
1315f2c45f1SShri Abhyankar   for (i = network->pStart; i < network->pEnd; i++) {
1325f2c45f1SShri Abhyankar     network->header[i].ndata = 0;
1335f2c45f1SShri Abhyankar     ndata = network->header[i].ndata;
1345f2c45f1SShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
1355f2c45f1SShri Abhyankar     network->header[i].offset[ndata] = 0;
1365f2c45f1SShri Abhyankar   }
137854ce69bSBarry Smith   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
1385f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
1395f2c45f1SShri Abhyankar }
1405f2c45f1SShri Abhyankar 
1415f2c45f1SShri Abhyankar #undef __FUNCT__
1425f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkRegisterComponent"
1435f2c45f1SShri Abhyankar /*@
1445f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
1455f2c45f1SShri Abhyankar 
1465f2c45f1SShri Abhyankar   Logically collective on DM
1475f2c45f1SShri Abhyankar 
1485f2c45f1SShri Abhyankar   Input Parameters
1495f2c45f1SShri Abhyankar + dm   - the network object
1505f2c45f1SShri Abhyankar . name - the component name
1515f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
1525f2c45f1SShri Abhyankar 
1535f2c45f1SShri Abhyankar    Output Parameters
1545f2c45f1SShri Abhyankar .   key - an integer key that defines the component
1555f2c45f1SShri Abhyankar 
1565f2c45f1SShri Abhyankar    Notes
1575f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
1585f2c45f1SShri Abhyankar 
1595f2c45f1SShri Abhyankar    Level: intermediate
1605f2c45f1SShri Abhyankar 
1615f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
1625f2c45f1SShri Abhyankar @*/
1635f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
1645f2c45f1SShri Abhyankar {
1655f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
1665f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
1675f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
1685f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
1695f2c45f1SShri Abhyankar   PetscInt              i;
1705f2c45f1SShri Abhyankar 
1715f2c45f1SShri Abhyankar   PetscFunctionBegin;
1725f2c45f1SShri Abhyankar 
1735f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
1745f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
1755f2c45f1SShri Abhyankar     if (flg) {
1765f2c45f1SShri Abhyankar       *key = i;
1775f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
1785f2c45f1SShri Abhyankar     }
1795f2c45f1SShri Abhyankar   }
1805f2c45f1SShri Abhyankar 
1815f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
1825f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
1835f2c45f1SShri Abhyankar   *key = network->ncomponent;
1845f2c45f1SShri Abhyankar   network->ncomponent++;
1855f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
1865f2c45f1SShri Abhyankar }
1875f2c45f1SShri Abhyankar 
1885f2c45f1SShri Abhyankar #undef __FUNCT__
1895f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVertexRange"
1905f2c45f1SShri Abhyankar /*@
1915f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
1925f2c45f1SShri Abhyankar 
1935f2c45f1SShri Abhyankar   Not Collective
1945f2c45f1SShri Abhyankar 
1955f2c45f1SShri Abhyankar   Input Parameters:
1965f2c45f1SShri Abhyankar + dm - The DMNetwork object
1975f2c45f1SShri Abhyankar 
1985f2c45f1SShri Abhyankar   Output Paramters:
1995f2c45f1SShri Abhyankar + vStart - The first vertex point
2005f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
2015f2c45f1SShri Abhyankar 
2025f2c45f1SShri Abhyankar   Level: intermediate
2035f2c45f1SShri Abhyankar 
2045f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
2055f2c45f1SShri Abhyankar @*/
2065f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
2075f2c45f1SShri Abhyankar {
2085f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
2095f2c45f1SShri Abhyankar 
2105f2c45f1SShri Abhyankar   PetscFunctionBegin;
2115f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
2125f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
2135f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2145f2c45f1SShri Abhyankar }
2155f2c45f1SShri Abhyankar 
2165f2c45f1SShri Abhyankar #undef __FUNCT__
2175f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetEdgeRange"
2185f2c45f1SShri Abhyankar /*@
2195f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
2205f2c45f1SShri Abhyankar 
2215f2c45f1SShri Abhyankar   Not Collective
2225f2c45f1SShri Abhyankar 
2235f2c45f1SShri Abhyankar   Input Parameters:
2245f2c45f1SShri Abhyankar + dm - The DMNetwork object
2255f2c45f1SShri Abhyankar 
2265f2c45f1SShri Abhyankar   Output Paramters:
2275f2c45f1SShri Abhyankar + eStart - The first edge point
2285f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
2295f2c45f1SShri Abhyankar 
2305f2c45f1SShri Abhyankar   Level: intermediate
2315f2c45f1SShri Abhyankar 
2325f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
2335f2c45f1SShri Abhyankar @*/
2345f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
2355f2c45f1SShri Abhyankar {
2365f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
2375f2c45f1SShri Abhyankar 
2385f2c45f1SShri Abhyankar   PetscFunctionBegin;
2395f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
2405f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
2415f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2425f2c45f1SShri Abhyankar }
2435f2c45f1SShri Abhyankar 
2445f2c45f1SShri Abhyankar #undef __FUNCT__
2455f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddComponent"
2465f2c45f1SShri Abhyankar /*@
2475f2c45f1SShri Abhyankar   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
2485f2c45f1SShri Abhyankar 
2495f2c45f1SShri Abhyankar   Not Collective
2505f2c45f1SShri Abhyankar 
2515f2c45f1SShri Abhyankar   Input Parameters:
2525f2c45f1SShri Abhyankar + dm           - The DMNetwork object
2535f2c45f1SShri Abhyankar . p            - vertex/edge point
2545f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
2555f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
2565f2c45f1SShri Abhyankar 
2575f2c45f1SShri Abhyankar   Level: intermediate
2585f2c45f1SShri Abhyankar 
2595f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
2605f2c45f1SShri Abhyankar @*/
2615f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
2625f2c45f1SShri Abhyankar {
2635f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
26443a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
2655f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
2665f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
2675f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
2685f2c45f1SShri Abhyankar 
2695f2c45f1SShri Abhyankar   PetscFunctionBegin;
27043a39a44SBarry Smith   header->size[header->ndata] = component->size;
27143a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
2725f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
2735f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
2745f2c45f1SShri Abhyankar 
2755f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
2765f2c45f1SShri Abhyankar   header->ndata++;
2775f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2785f2c45f1SShri Abhyankar }
2795f2c45f1SShri Abhyankar 
2805f2c45f1SShri Abhyankar #undef __FUNCT__
2815f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetNumComponents"
2825f2c45f1SShri Abhyankar /*@
2835f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
2845f2c45f1SShri Abhyankar 
2855f2c45f1SShri Abhyankar   Not Collective
2865f2c45f1SShri Abhyankar 
2875f2c45f1SShri Abhyankar   Input Parameters:
2885f2c45f1SShri Abhyankar + dm - The DMNetwork object
2895f2c45f1SShri Abhyankar . p  - vertex/edge point
2905f2c45f1SShri Abhyankar 
2915f2c45f1SShri Abhyankar   Output Parameters:
2925f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
2935f2c45f1SShri Abhyankar 
2945f2c45f1SShri Abhyankar   Level: intermediate
2955f2c45f1SShri Abhyankar 
2965f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
2975f2c45f1SShri Abhyankar @*/
2985f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
2995f2c45f1SShri Abhyankar {
3005f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3015f2c45f1SShri Abhyankar   PetscInt       offset;
3025f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3035f2c45f1SShri Abhyankar 
3045f2c45f1SShri Abhyankar   PetscFunctionBegin;
3055f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
3065f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
3075f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3085f2c45f1SShri Abhyankar }
3095f2c45f1SShri Abhyankar 
3105f2c45f1SShri Abhyankar #undef __FUNCT__
3115f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetComponentTypeOffset"
3125f2c45f1SShri Abhyankar /*@
3135f2c45f1SShri Abhyankar   DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the
3145f2c45f1SShri Abhyankar                                     component value from the component data array
3155f2c45f1SShri Abhyankar 
3165f2c45f1SShri Abhyankar   Not Collective
3175f2c45f1SShri Abhyankar 
3185f2c45f1SShri Abhyankar   Input Parameters:
3195f2c45f1SShri Abhyankar + dm      - The DMNetwork object
3205f2c45f1SShri Abhyankar . p       - vertex/edge point
3215f2c45f1SShri Abhyankar - compnum - component number
3225f2c45f1SShri Abhyankar 
3235f2c45f1SShri Abhyankar   Output Parameters:
3245f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component
3255f2c45f1SShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
3265f2c45f1SShri Abhyankar 
3275f2c45f1SShri Abhyankar   Notes:
3285f2c45f1SShri Abhyankar   Typical usage:
3295f2c45f1SShri Abhyankar 
3305f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
3315f2c45f1SShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
3325f2c45f1SShri Abhyankar   Loop over vertices or edges
3335f2c45f1SShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
3345f2c45f1SShri Abhyankar     Loop over numcomps
3355f2c45f1SShri Abhyankar       DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset);
3365f2c45f1SShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
3375f2c45f1SShri Abhyankar 
3385f2c45f1SShri Abhyankar   Level: intermediate
3395f2c45f1SShri Abhyankar 
3405f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
3415f2c45f1SShri Abhyankar @*/
3425f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
3435f2c45f1SShri Abhyankar {
3445f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
3455f2c45f1SShri Abhyankar   PetscInt                 offsetp;
3465f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
3475f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
3485f2c45f1SShri Abhyankar 
3495f2c45f1SShri Abhyankar   PetscFunctionBegin;
3505f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
3515f2c45f1SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
3525f2c45f1SShri Abhyankar   *compkey = header->key[compnum];
3535f2c45f1SShri Abhyankar   *offset  = offsetp+network->dataheadersize+header->offset[compnum];
3545f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3555f2c45f1SShri Abhyankar }
3565f2c45f1SShri Abhyankar 
3575f2c45f1SShri Abhyankar #undef __FUNCT__
3585f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVariableOffset"
3595f2c45f1SShri Abhyankar /*@
3605f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
3615f2c45f1SShri Abhyankar 
3625f2c45f1SShri Abhyankar   Not Collective
3635f2c45f1SShri Abhyankar 
3645f2c45f1SShri Abhyankar   Input Parameters:
3655f2c45f1SShri Abhyankar + dm     - The DMNetwork object
3665f2c45f1SShri Abhyankar - p      - the edge/vertex point
3675f2c45f1SShri Abhyankar 
3685f2c45f1SShri Abhyankar   Output Parameters:
3695f2c45f1SShri Abhyankar . offset - the offset
3705f2c45f1SShri Abhyankar 
3715f2c45f1SShri Abhyankar   Level: intermediate
3725f2c45f1SShri Abhyankar 
3735f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
3745f2c45f1SShri Abhyankar @*/
3755f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
3765f2c45f1SShri Abhyankar {
3775f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3785f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3795f2c45f1SShri Abhyankar 
3805f2c45f1SShri Abhyankar   PetscFunctionBegin;
3815f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DofSection,p,offset);CHKERRQ(ierr);
3825f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3835f2c45f1SShri Abhyankar }
3845f2c45f1SShri Abhyankar 
3855f2c45f1SShri Abhyankar #undef __FUNCT__
3865f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVariableGlobalOffset"
3875f2c45f1SShri Abhyankar /*@
3885f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
3895f2c45f1SShri Abhyankar 
3905f2c45f1SShri Abhyankar   Not Collective
3915f2c45f1SShri Abhyankar 
3925f2c45f1SShri Abhyankar   Input Parameters:
3935f2c45f1SShri Abhyankar + dm      - The DMNetwork object
3945f2c45f1SShri Abhyankar - p       - the edge/vertex point
3955f2c45f1SShri Abhyankar 
3965f2c45f1SShri Abhyankar   Output Parameters:
3975f2c45f1SShri Abhyankar . offsetg - the offset
3985f2c45f1SShri Abhyankar 
3995f2c45f1SShri Abhyankar   Level: intermediate
4005f2c45f1SShri Abhyankar 
4015f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
4025f2c45f1SShri Abhyankar @*/
4035f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
4045f2c45f1SShri Abhyankar {
4055f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4065f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4075f2c45f1SShri Abhyankar 
4085f2c45f1SShri Abhyankar   PetscFunctionBegin;
4095f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);CHKERRQ(ierr);
4105f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4115f2c45f1SShri Abhyankar }
4125f2c45f1SShri Abhyankar 
4135f2c45f1SShri Abhyankar #undef __FUNCT__
4145f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddNumVariables"
4155f2c45f1SShri Abhyankar /*@
4165f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
4175f2c45f1SShri Abhyankar 
4185f2c45f1SShri Abhyankar   Not Collective
4195f2c45f1SShri Abhyankar 
4205f2c45f1SShri Abhyankar   Input Parameters:
4215f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
4225f2c45f1SShri Abhyankar . p    - the vertex/edge point
4235f2c45f1SShri Abhyankar - nvar - number of additional variables
4245f2c45f1SShri Abhyankar 
4255f2c45f1SShri Abhyankar   Level: intermediate
4265f2c45f1SShri Abhyankar 
4275f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
4285f2c45f1SShri Abhyankar @*/
4295f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
4305f2c45f1SShri Abhyankar {
4315f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4325f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4335f2c45f1SShri Abhyankar 
4345f2c45f1SShri Abhyankar   PetscFunctionBegin;
4355f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
4365f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4375f2c45f1SShri Abhyankar }
4385f2c45f1SShri Abhyankar 
4395f2c45f1SShri Abhyankar #undef __FUNCT__
44027f51fceSHong Zhang #define __FUNCT__ "DMNetworkGetNumVariables"
44127f51fceSHong Zhang /*@
44227f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
44327f51fceSHong Zhang 
44427f51fceSHong Zhang   Not Collective
44527f51fceSHong Zhang 
44627f51fceSHong Zhang   Input Parameters:
44727f51fceSHong Zhang + dm   - The DMNetworkObject
44827f51fceSHong Zhang - p    - the vertex/edge point
44927f51fceSHong Zhang 
45027f51fceSHong Zhang   Output Parameters:
45127f51fceSHong Zhang . nvar - number of variables
45227f51fceSHong Zhang 
45327f51fceSHong Zhang   Level: intermediate
45427f51fceSHong Zhang 
45527f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
45627f51fceSHong Zhang @*/
45727f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
45827f51fceSHong Zhang {
45927f51fceSHong Zhang   PetscErrorCode ierr;
46027f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
46127f51fceSHong Zhang 
46227f51fceSHong Zhang   PetscFunctionBegin;
46327f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
46427f51fceSHong Zhang   PetscFunctionReturn(0);
46527f51fceSHong Zhang }
46627f51fceSHong Zhang 
46727f51fceSHong Zhang #undef __FUNCT__
468662b5b01SHong Zhang #define __FUNCT__ "DMNetworkSetNumVariables"
4695f2c45f1SShri Abhyankar /*@
4705f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
4715f2c45f1SShri Abhyankar 
4725f2c45f1SShri Abhyankar   Not Collective
4735f2c45f1SShri Abhyankar 
4745f2c45f1SShri Abhyankar   Input Parameters:
4755f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
4765f2c45f1SShri Abhyankar . p    - the vertex/edge point
4775f2c45f1SShri Abhyankar - nvar - number of variables
4785f2c45f1SShri Abhyankar 
4795f2c45f1SShri Abhyankar   Level: intermediate
4805f2c45f1SShri Abhyankar 
4815f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
4825f2c45f1SShri Abhyankar @*/
4835f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
4845f2c45f1SShri Abhyankar {
4855f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4865f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4875f2c45f1SShri Abhyankar 
4885f2c45f1SShri Abhyankar   PetscFunctionBegin;
4895f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
4905f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4915f2c45f1SShri Abhyankar }
4925f2c45f1SShri Abhyankar 
4935f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
4945f2c45f1SShri Abhyankar    function is called during DMSetUp() */
4955f2c45f1SShri Abhyankar #undef __FUNCT__
4965f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkComponentSetUp"
4975f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
4985f2c45f1SShri Abhyankar {
4995f2c45f1SShri Abhyankar   PetscErrorCode              ierr;
5005f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5015f2c45f1SShri Abhyankar   PetscInt                    arr_size;
5025f2c45f1SShri Abhyankar   PetscInt                    p,offset,offsetp;
5035f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
5045f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
5055f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType      *componentdataarray;
5065f2c45f1SShri Abhyankar   PetscInt ncomp, i;
5075f2c45f1SShri Abhyankar 
5085f2c45f1SShri Abhyankar   PetscFunctionBegin;
5095f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
5105f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
51175b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
5125f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
5135f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
5145f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
5155f2c45f1SShri Abhyankar     /* Copy header */
5165f2c45f1SShri Abhyankar     header = &network->header[p];
517302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
5185f2c45f1SShri Abhyankar     /* Copy data */
5195f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
5205f2c45f1SShri Abhyankar     ncomp = header->ndata;
5215f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
5225f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
523302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
5245f2c45f1SShri Abhyankar     }
5255f2c45f1SShri Abhyankar   }
5265f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5275f2c45f1SShri Abhyankar }
5285f2c45f1SShri Abhyankar 
5295f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
5305f2c45f1SShri Abhyankar #undef __FUNCT__
5315f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkVariablesSetUp"
5325f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
5335f2c45f1SShri Abhyankar {
5345f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5355f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5365f2c45f1SShri Abhyankar 
5375f2c45f1SShri Abhyankar   PetscFunctionBegin;
5385f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
5395f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5405f2c45f1SShri Abhyankar }
5415f2c45f1SShri Abhyankar 
5425f2c45f1SShri Abhyankar #undef __FUNCT__
5435f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetComponentDataArray"
5445f2c45f1SShri Abhyankar /*@C
5455f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
5465f2c45f1SShri Abhyankar 
5475f2c45f1SShri Abhyankar   Not Collective
5485f2c45f1SShri Abhyankar 
5495f2c45f1SShri Abhyankar   Input Parameters:
5505f2c45f1SShri Abhyankar . dm - The DMNetwork Object
5515f2c45f1SShri Abhyankar 
5525f2c45f1SShri Abhyankar   Output Parameters:
5535f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
5545f2c45f1SShri Abhyankar 
5555f2c45f1SShri Abhyankar   Level: intermediate
5565f2c45f1SShri Abhyankar 
5575f2c45f1SShri Abhyankar .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents
5585f2c45f1SShri Abhyankar @*/
5595f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
5605f2c45f1SShri Abhyankar {
5615f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5625f2c45f1SShri Abhyankar 
5635f2c45f1SShri Abhyankar   PetscFunctionBegin;
5645f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
5655f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5665f2c45f1SShri Abhyankar }
5675f2c45f1SShri Abhyankar 
5685f2c45f1SShri Abhyankar #undef __FUNCT__
5695f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkDistribute"
5705f2c45f1SShri Abhyankar /*@
5715f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
5725f2c45f1SShri Abhyankar 
5735f2c45f1SShri Abhyankar   Collective
5745f2c45f1SShri Abhyankar 
5755f2c45f1SShri Abhyankar   Input Parameter:
5765f2c45f1SShri Abhyankar + oldDM - the original DMNetwork object
5775f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
5785f2c45f1SShri Abhyankar 
5795f2c45f1SShri Abhyankar   Output Parameter:
5805f2c45f1SShri Abhyankar . distDM - the distributed DMNetwork object
5815f2c45f1SShri Abhyankar 
5825f2c45f1SShri Abhyankar   Notes:
5835f2c45f1SShri Abhyankar   This routine should be called only when using multiple processors.
5845f2c45f1SShri Abhyankar 
5858b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
5865f2c45f1SShri Abhyankar 
5875f2c45f1SShri Abhyankar   Level: intermediate
5885f2c45f1SShri Abhyankar 
5895f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
5905f2c45f1SShri Abhyankar @*/
59180cf41d5SMatthew G. Knepley PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM)
5925f2c45f1SShri Abhyankar {
5935f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5945f2c45f1SShri Abhyankar   DM_Network     *oldDMnetwork = (DM_Network*)oldDM->data;
5955f2c45f1SShri Abhyankar   PetscSF        pointsf;
5965f2c45f1SShri Abhyankar   DM             newDM;
5975f2c45f1SShri Abhyankar   DM_Network     *newDMnetwork;
5985f2c45f1SShri Abhyankar 
5995f2c45f1SShri Abhyankar   PetscFunctionBegin;
6005f2c45f1SShri Abhyankar   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr);
6015f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
6025f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
6035f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
60480cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
6055f2c45f1SShri Abhyankar   /* Distribute dof section */
6065f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr);
6075f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
6085f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr);
6095f2c45f1SShri Abhyankar   /* Distribute data and associated section */
61031da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
6115f2c45f1SShri Abhyankar   /* Destroy point SF */
6125f2c45f1SShri Abhyankar   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
6135f2c45f1SShri Abhyankar 
6145f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
6155f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
6165f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
6175f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
6185f2c45f1SShri Abhyankar   newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart;
6195f2c45f1SShri Abhyankar   newDMnetwork->NNodes = oldDMnetwork->NNodes;
6205f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
6215f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
6225f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
6235f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
6245f2c45f1SShri Abhyankar 
6255f2c45f1SShri Abhyankar   *distDM = newDM;
6265f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6275f2c45f1SShri Abhyankar }
6285f2c45f1SShri Abhyankar 
6295f2c45f1SShri Abhyankar #undef __FUNCT__
6305f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetSupportingEdges"
6315f2c45f1SShri Abhyankar /*@C
6325f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
6335f2c45f1SShri Abhyankar 
6345f2c45f1SShri Abhyankar   Not Collective
6355f2c45f1SShri Abhyankar 
6365f2c45f1SShri Abhyankar   Input Parameters:
6375f2c45f1SShri Abhyankar + dm - The DMNetwork object
6385f2c45f1SShri Abhyankar - p  - the vertex point
6395f2c45f1SShri Abhyankar 
6405f2c45f1SShri Abhyankar   Output Paramters:
6415f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
6425f2c45f1SShri Abhyankar - edges  - List of edge points
6435f2c45f1SShri Abhyankar 
6445f2c45f1SShri Abhyankar   Level: intermediate
6455f2c45f1SShri Abhyankar 
6465f2c45f1SShri Abhyankar   Fortran Notes:
6475f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
6485f2c45f1SShri Abhyankar   include petsc.h90 in your code.
6495f2c45f1SShri Abhyankar 
6505f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
6515f2c45f1SShri Abhyankar @*/
6525f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
6535f2c45f1SShri Abhyankar {
6545f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6555f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6565f2c45f1SShri Abhyankar 
6575f2c45f1SShri Abhyankar   PetscFunctionBegin;
6585f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
6595f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
6605f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6615f2c45f1SShri Abhyankar }
6625f2c45f1SShri Abhyankar 
6635f2c45f1SShri Abhyankar #undef __FUNCT__
6645f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetConnectedNodes"
6655f2c45f1SShri Abhyankar /*@C
666596e729fSHong Zhang   DMNetworkGetConnectedNodes - Return the connected vertices for this edge point
6675f2c45f1SShri Abhyankar 
6685f2c45f1SShri Abhyankar   Not Collective
6695f2c45f1SShri Abhyankar 
6705f2c45f1SShri Abhyankar   Input Parameters:
6715f2c45f1SShri Abhyankar + dm - The DMNetwork object
6725f2c45f1SShri Abhyankar - p  - the edge point
6735f2c45f1SShri Abhyankar 
6745f2c45f1SShri Abhyankar   Output Paramters:
6755f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
6765f2c45f1SShri Abhyankar 
6775f2c45f1SShri Abhyankar   Level: intermediate
6785f2c45f1SShri Abhyankar 
6795f2c45f1SShri Abhyankar   Fortran Notes:
6805f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
6815f2c45f1SShri Abhyankar   include petsc.h90 in your code.
6825f2c45f1SShri Abhyankar 
6835f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
6845f2c45f1SShri Abhyankar @*/
6855f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
6865f2c45f1SShri Abhyankar {
6875f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6885f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6895f2c45f1SShri Abhyankar 
6905f2c45f1SShri Abhyankar   PetscFunctionBegin;
6915f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
6925f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6935f2c45f1SShri Abhyankar }
6945f2c45f1SShri Abhyankar 
6955f2c45f1SShri Abhyankar #undef __FUNCT__
6965f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkIsGhostVertex"
6975f2c45f1SShri Abhyankar /*@
6985f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
6995f2c45f1SShri Abhyankar 
7005f2c45f1SShri Abhyankar   Not Collective
7015f2c45f1SShri Abhyankar 
7025f2c45f1SShri Abhyankar   Input Parameters:
7035f2c45f1SShri Abhyankar + dm - The DMNetwork object
7045f2c45f1SShri Abhyankar . p  - the vertex point
7055f2c45f1SShri Abhyankar 
7065f2c45f1SShri Abhyankar   Output Parameter:
7075f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
7085f2c45f1SShri Abhyankar 
7095f2c45f1SShri Abhyankar   Level: intermediate
7105f2c45f1SShri Abhyankar 
7115f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
7125f2c45f1SShri Abhyankar @*/
7135f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
7145f2c45f1SShri Abhyankar {
7155f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7165f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7175f2c45f1SShri Abhyankar   PetscInt       offsetg;
7185f2c45f1SShri Abhyankar   PetscSection   sectiong;
7195f2c45f1SShri Abhyankar 
7205f2c45f1SShri Abhyankar   PetscFunctionBegin;
7215f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
7225f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
7235f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
7245f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
7255f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7265f2c45f1SShri Abhyankar }
7275f2c45f1SShri Abhyankar 
7285f2c45f1SShri Abhyankar #undef __FUNCT__
7295f2c45f1SShri Abhyankar #define __FUNCT__ "DMSetUp_Network"
7305f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
7315f2c45f1SShri Abhyankar {
7325f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7335f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
7345f2c45f1SShri Abhyankar 
7355f2c45f1SShri Abhyankar   PetscFunctionBegin;
7365f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
7375f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
7385f2c45f1SShri Abhyankar 
7395f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
7405f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
7415f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7425f2c45f1SShri Abhyankar }
7435f2c45f1SShri Abhyankar 
7445f2c45f1SShri Abhyankar #undef __FUNCT__
74517df6e9eSHong Zhang #define __FUNCT__ "DMNetworkHasJacobian"
7461ad426b7SHong Zhang /*@
74717df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
7481ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
7491ad426b7SHong Zhang 
7501ad426b7SHong Zhang     Collective
7511ad426b7SHong Zhang 
7521ad426b7SHong Zhang     Input Parameters:
753*83b2e829SHong Zhang +   dm - The DMNetwork object
754*83b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
755*83b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
7561ad426b7SHong Zhang 
7571ad426b7SHong Zhang     Level: intermediate
7581ad426b7SHong Zhang 
7591ad426b7SHong Zhang @*/
760*83b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
7611ad426b7SHong Zhang {
7621ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
7631ad426b7SHong Zhang 
7641ad426b7SHong Zhang   PetscFunctionBegin;
765*83b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
766*83b2e829SHong Zhang   network->userVertexJacobian = vflg;
7671ad426b7SHong Zhang   PetscFunctionReturn(0);
7681ad426b7SHong Zhang }
7691ad426b7SHong Zhang 
7701ad426b7SHong Zhang #undef __FUNCT__
771*83b2e829SHong Zhang #define __FUNCT__ "DMNetworkEdgeSetMatrix"
7721ad426b7SHong Zhang /*@
773*83b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
774*83b2e829SHong Zhang 
775*83b2e829SHong Zhang     Not Collective
776*83b2e829SHong Zhang 
777*83b2e829SHong Zhang     Input Parameters:
778*83b2e829SHong Zhang +   dm - The DMNetwork object
779*83b2e829SHong Zhang .   p  - the edge point
780*83b2e829SHong Zhang -   J - array of Jacobian submatrices:
781*83b2e829SHong Zhang         J[0]: this edge;
782*83b2e829SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes();
783*83b2e829SHong Zhang 
784*83b2e829SHong Zhang     Level: intermediate
785*83b2e829SHong Zhang 
786*83b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
787*83b2e829SHong Zhang @*/
788*83b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
789*83b2e829SHong Zhang {
790*83b2e829SHong Zhang   PetscErrorCode ierr;
791*83b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
792*83b2e829SHong Zhang 
793*83b2e829SHong Zhang   PetscFunctionBegin;
794*83b2e829SHong Zhang   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkCreateJacobian() collectively before calling DMNetworkElementSetMatrix");
795*83b2e829SHong Zhang   if (!network->Je) {
796*83b2e829SHong Zhang     ierr = PetscCalloc1(network->nEdges,&network->Je);CHKERRQ(ierr);
797*83b2e829SHong Zhang   }
798*83b2e829SHong Zhang   network->Je[p] = J[0];
799*83b2e829SHong Zhang   PetscFunctionReturn(0);
800*83b2e829SHong Zhang }
801*83b2e829SHong Zhang 
802*83b2e829SHong Zhang #undef __FUNCT__
803*83b2e829SHong Zhang #define __FUNCT__ "DMNetworkVetexSetMatrix"
804*83b2e829SHong Zhang /*@
805*83b2e829SHong Zhang     DMNetworkVetexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
8061ad426b7SHong Zhang 
8071ad426b7SHong Zhang     Not Collective
8081ad426b7SHong Zhang 
8091ad426b7SHong Zhang     Input Parameters:
8101ad426b7SHong Zhang +   dm - The DMNetwork object
8111ad426b7SHong Zhang .   p  - the vertex point
812*83b2e829SHong Zhang -   J - array of Jacobian submatrices:
8131ad426b7SHong Zhang 
8141ad426b7SHong Zhang     Level: intermediate
8151ad426b7SHong Zhang 
816*83b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
8171ad426b7SHong Zhang @*/
818*83b2e829SHong Zhang PetscErrorCode DMNetworkVetexSetMatrix(DM dm,PetscInt p,Mat J[])
8195f2c45f1SShri Abhyankar {
8205f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8215f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
8225f2c45f1SShri Abhyankar 
8235f2c45f1SShri Abhyankar   PetscFunctionBegin;
824*83b2e829SHong Zhang   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkCreateJacobian() collectively before calling DMNetworkElementSetMatrix");
825*83b2e829SHong Zhang   if (!network->Jv) {
826*83b2e829SHong Zhang     ierr = PetscCalloc1(network->nNodes,&network->Jv);CHKERRQ(ierr);
8271ad426b7SHong Zhang   }
828*83b2e829SHong Zhang   network->Jv[p] = J[0];
8291ad426b7SHong Zhang   PetscFunctionReturn(0);
8301ad426b7SHong Zhang }
8311ad426b7SHong Zhang 
8321ad426b7SHong Zhang #undef __FUNCT__
8331ad426b7SHong Zhang #define __FUNCT__ "DMCreateMatrix_Network"
8341ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
8351ad426b7SHong Zhang {
8361ad426b7SHong Zhang   PetscErrorCode ierr;
8371ad426b7SHong Zhang   DM_Network     *network = (DM_Network*) dm->data;
838bfbc38dcSHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,rend,row,row_e,nrows,localSize;
839bfbc38dcSHong Zhang   PetscInt       cstart,ncols,col,j,e,v,*dnz,*onz,*dnzu,*onzu;
8401ad426b7SHong Zhang   const PetscInt *cols;
84117df6e9eSHong Zhang   PetscScalar    *zeros,zero = 0.0;
8421ad426b7SHong Zhang   PetscBool      ghost;
84317df6e9eSHong Zhang   Mat            Je;
844bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
845bfbc38dcSHong Zhang   PetscInt       nedges;
846bfbc38dcSHong Zhang   const PetscInt *edges;
8471ad426b7SHong Zhang 
8481ad426b7SHong Zhang   PetscFunctionBegin;
849*83b2e829SHong Zhang   if (!network->userEdgeJacobian && !network->userVertexJacobian) { /* user does not provide Jacobian blocks */
850bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
851bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
8521ad426b7SHong Zhang     PetscFunctionReturn(0);
8531ad426b7SHong Zhang   }
8541ad426b7SHong Zhang 
855bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
8562a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
857bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
858bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
8592a945128SHong Zhang 
8602a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
8612a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
8621ad426b7SHong Zhang 
863bfbc38dcSHong Zhang   /* Preallocation - submatrix for an element (edge/vertex) is allocated as a dense block, see DMCreateMatrix_Plex() */
864bfbc38dcSHong Zhang   ierr = PetscCalloc4(localSize,&dnz,localSize,&onz,localSize,&dnzu,localSize,&onzu);CHKERRQ(ierr);
865bfbc38dcSHong Zhang   ierr = DMPlexPreallocateOperator(network->plex,1,dnz,onz,dnzu,onzu,*J,PETSC_FALSE);CHKERRQ(ierr);
866bfbc38dcSHong Zhang   ierr = PetscFree4(dnz,onz,dnzu,onzu);CHKERRQ(ierr);
867bfbc38dcSHong Zhang 
868bfbc38dcSHong Zhang   /* Set matrix entries for edges */
869*83b2e829SHong Zhang   /*------------------------------*/
870bfbc38dcSHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
8711ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
872bfbc38dcSHong Zhang     /* Get row indices */
8731ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
87417df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
8751ad426b7SHong Zhang 
87617df6e9eSHong Zhang     PetscInt    rows[nrows],*cols_tmp;
87717df6e9eSHong Zhang     for (j=0; j<nrows; j++) rows[j] = j + rstart;
8781ad426b7SHong Zhang 
879bfbc38dcSHong Zhang     /* Set matrix entries for conntected vertices */
88017df6e9eSHong Zhang     const PetscInt    *cone;
88117df6e9eSHong Zhang     ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
88217df6e9eSHong Zhang 
883bfbc38dcSHong Zhang     for (v=0; v<2; v++) {
884bfbc38dcSHong Zhang       ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
885bfbc38dcSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
886bfbc38dcSHong Zhang       if (ghost) cstart = -(cstart + 1); /* Convert to actual global offset for ghost nodes */
887bfbc38dcSHong Zhang       ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
88817df6e9eSHong Zhang       ierr = PetscCalloc2(ncols,&cols_tmp,nrows*ncols,&zeros);CHKERRQ(ierr);
889bfbc38dcSHong Zhang       for (j=0; j<ncols; j++) cols_tmp[j] = j+ cstart;
89017df6e9eSHong Zhang 
89117df6e9eSHong Zhang       ierr = MatSetValues(*J,nrows,rows,ncols,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr);
89217df6e9eSHong Zhang       ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr);
893bfbc38dcSHong Zhang     }
89417df6e9eSHong Zhang 
895bfbc38dcSHong Zhang     /* Set matrix entries for edge self */
896*83b2e829SHong Zhang     Je = network->Je[e];
897*83b2e829SHong Zhang     if (Je) { /* use user-provided Je */
898*83b2e829SHong Zhang       PetscInt M,N;
899*83b2e829SHong Zhang       ierr = MatGetSize(Je,&M,&N);CHKERRQ(ierr);
900*83b2e829SHong Zhang       if (nrows != M || nrows != N) SETERRQ3(PetscObjectComm((PetscObject)Je),PETSC_ERR_USER,"%D must equal %D and %D",rend-rstart,M,N);
901*83b2e829SHong Zhang 
902bfbc38dcSHong Zhang       rend = rstart + nrows;
9031ad426b7SHong Zhang       for (row=rstart; row<rend; row++) {
9041ad426b7SHong Zhang         row_e = row - rstart;
90517df6e9eSHong Zhang         ierr = MatGetRow(Je,row_e,&ncols,&cols,NULL);CHKERRQ(ierr);
9061ad426b7SHong Zhang         for (j=0; j<ncols; j++) {
90717df6e9eSHong Zhang           col = cols[j] + rstart;
90817df6e9eSHong Zhang           ierr = MatSetValues(*J,1,&row,1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
9091ad426b7SHong Zhang         }
91017df6e9eSHong Zhang         ierr = MatRestoreRow(Je,row_e,&ncols,&cols,NULL);CHKERRQ(ierr);
9111ad426b7SHong Zhang       }
912*83b2e829SHong Zhang     } else { /* set a dense block */
913*83b2e829SHong Zhang       PetscInt *cols_tmp,cstart;
914*83b2e829SHong Zhang       ierr = PetscCalloc2(nrows,&cols_tmp,nrows*nrows,&zeros);CHKERRQ(ierr);
915*83b2e829SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,e,&cstart);CHKERRQ(ierr);
916*83b2e829SHong Zhang       for (j=0; j<nrows; j++) cols_tmp[j] = j+ cstart;
917*83b2e829SHong Zhang       ierr = MatSetValues(*J,nrows,rows,nrows,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr);
918*83b2e829SHong Zhang       ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr);
919*83b2e829SHong Zhang     }
9201ad426b7SHong Zhang   }
9211ad426b7SHong Zhang 
922bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
923*83b2e829SHong Zhang   /*---------------------------------*/
924bfbc38dcSHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
9251ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
926bfbc38dcSHong Zhang     /* Get row indices */
9272a945128SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
928596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
9292a945128SHong Zhang     if (ghost) rstart = -(rstart + 1); /* Convert to actual global offset for ghost nodes */
930596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
931596e729fSHong Zhang 
932596e729fSHong Zhang     PetscInt    rows[nrows];
933596e729fSHong Zhang     for (j=0; j<nrows; j++) rows[j] = j + rstart;
934596e729fSHong Zhang 
935bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
936596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
937596e729fSHong Zhang 
938596e729fSHong Zhang     for (e=0; e<nedges; e++) {
939bfbc38dcSHong Zhang       /* Supporting edges */
940bfbc38dcSHong Zhang       PetscInt *cols_tmp;
941596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
942596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
943596e729fSHong Zhang 
944596e729fSHong Zhang       ierr = PetscCalloc2(ncols,&cols_tmp,nrows*ncols,&zeros);CHKERRQ(ierr);
945596e729fSHong Zhang       for (j=0; j<ncols; j++) cols_tmp[j] = j+ cstart;
946596e729fSHong Zhang       ierr = MatSetValues(*J,nrows,rows,ncols,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr);
947596e729fSHong Zhang       ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr);
948596e729fSHong Zhang 
949bfbc38dcSHong Zhang       /* Connected vertices */
95044aca652SHong Zhang       const PetscInt *cone;
95144aca652SHong Zhang       PetscInt       vc;
95244aca652SHong Zhang       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
9532a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
9542a945128SHong Zhang 
95544aca652SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost);CHKERRQ(ierr);
95644aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
95744aca652SHong Zhang       if (ghost) cstart = -(cstart + 1); /* Convert to actual global offset for ghost nodes */
95844aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
95944aca652SHong Zhang 
96044aca652SHong Zhang       ierr = PetscCalloc2(ncols,&cols_tmp,nrows*ncols,&zeros);CHKERRQ(ierr);
96144aca652SHong Zhang       for (j=0; j<ncols; j++) cols_tmp[j] = j+ cstart;
96244aca652SHong Zhang       ierr = MatSetValues(*J,nrows,rows,ncols,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr);
96344aca652SHong Zhang       ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr);
964596e729fSHong Zhang     }
965596e729fSHong Zhang 
966bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
9671ad426b7SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
9681ad426b7SHong Zhang     if (!ghost) {
969596e729fSHong Zhang       PetscInt *cols_tmp,cstart;
970596e729fSHong Zhang       ierr = PetscCalloc2(nrows,&cols_tmp,nrows*nrows,&zeros);CHKERRQ(ierr);
971596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
972596e729fSHong Zhang       for (j=0; j<nrows; j++) cols_tmp[j] = j+ cstart;
973596e729fSHong Zhang       ierr = MatSetValues(*J,nrows,rows,nrows,cols_tmp,zeros,INSERT_VALUES);CHKERRQ(ierr);
974596e729fSHong Zhang       ierr = PetscFree2(cols_tmp,zeros);CHKERRQ(ierr);
9751ad426b7SHong Zhang     }
9761ad426b7SHong Zhang   }
9771ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9781ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
9792a945128SHong Zhang #if 0
980bfbc38dcSHong Zhang   printf("\nMatrix J:\n");
9812a945128SHong Zhang   ierr = MatView(*J,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
9822a945128SHong Zhang #endif
9835f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
9845f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9855f2c45f1SShri Abhyankar }
9865f2c45f1SShri Abhyankar 
9875f2c45f1SShri Abhyankar #undef __FUNCT__
9885f2c45f1SShri Abhyankar #define __FUNCT__ "DMDestroy_Network"
9895f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
9905f2c45f1SShri Abhyankar {
9915f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9925f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
993*83b2e829SHong Zhang   PetscInt       i;
9945f2c45f1SShri Abhyankar 
9955f2c45f1SShri Abhyankar   PetscFunctionBegin;
9968415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
997*83b2e829SHong Zhang   if (network->Je) {
998*83b2e829SHong Zhang     for (i=0; i<network->nEdges; i++) {
999*83b2e829SHong Zhang       ierr = MatDestroy(&network->Je[i]);CHKERRQ(ierr);
10001ad426b7SHong Zhang     }
1001*83b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
1002*83b2e829SHong Zhang   }
1003*83b2e829SHong Zhang   if (network->Jv) {
1004*83b2e829SHong Zhang     for (i=0; i<network->nNodes; i++) {
1005*83b2e829SHong Zhang       ierr = MatDestroy(&network->Jv[i]);CHKERRQ(ierr);
1006*83b2e829SHong Zhang     }
1007*83b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
10081ad426b7SHong Zhang   }
10095f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
10105f2c45f1SShri Abhyankar   network->edges = NULL;
10115f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
10125f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
1013*83b2e829SHong Zhang 
10145f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
10155f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
10165f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
10175f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
10185f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10195f2c45f1SShri Abhyankar }
10205f2c45f1SShri Abhyankar 
10215f2c45f1SShri Abhyankar #undef __FUNCT__
10225f2c45f1SShri Abhyankar #define __FUNCT__ "DMView_Network"
10235f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
10245f2c45f1SShri Abhyankar {
10255f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10265f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
10275f2c45f1SShri Abhyankar 
10285f2c45f1SShri Abhyankar   PetscFunctionBegin;
10295f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
10305f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10315f2c45f1SShri Abhyankar }
10325f2c45f1SShri Abhyankar 
10335f2c45f1SShri Abhyankar #undef __FUNCT__
10345f2c45f1SShri Abhyankar #define __FUNCT__ "DMGlobalToLocalBegin_Network"
10355f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
10365f2c45f1SShri Abhyankar {
10375f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10385f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
10395f2c45f1SShri Abhyankar 
10405f2c45f1SShri Abhyankar   PetscFunctionBegin;
10415f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
10425f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10435f2c45f1SShri Abhyankar }
10445f2c45f1SShri Abhyankar 
10455f2c45f1SShri Abhyankar #undef __FUNCT__
10465f2c45f1SShri Abhyankar #define __FUNCT__ "DMGlobalToLocalEnd_Network"
10475f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
10485f2c45f1SShri Abhyankar {
10495f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10505f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
10515f2c45f1SShri Abhyankar 
10525f2c45f1SShri Abhyankar   PetscFunctionBegin;
10535f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
10545f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10555f2c45f1SShri Abhyankar }
10565f2c45f1SShri Abhyankar 
10575f2c45f1SShri Abhyankar #undef __FUNCT__
10585f2c45f1SShri Abhyankar #define __FUNCT__ "DMLocalToGlobalBegin_Network"
10595f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
10605f2c45f1SShri Abhyankar {
10615f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10625f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
10635f2c45f1SShri Abhyankar 
10645f2c45f1SShri Abhyankar   PetscFunctionBegin;
10655f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
10665f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10675f2c45f1SShri Abhyankar }
10685f2c45f1SShri Abhyankar 
10695f2c45f1SShri Abhyankar #undef __FUNCT__
10705f2c45f1SShri Abhyankar #define __FUNCT__ "DMLocalToGlobalEnd_Network"
10715f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
10725f2c45f1SShri Abhyankar {
10735f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10745f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
10755f2c45f1SShri Abhyankar 
10765f2c45f1SShri Abhyankar   PetscFunctionBegin;
10775f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
10785f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10795f2c45f1SShri Abhyankar }
1080