xref: /petsc/src/dm/impls/network/network.c (revision 840c2264659c88cc434f047d917da87fb6fa55c4)
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;
270fa58f0a9SHong Zhang   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);
271fa58f0a9SHong Zhang 
27243a39a44SBarry Smith   header->size[header->ndata] = component->size;
27343a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
2745f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
2755f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
2765f2c45f1SShri Abhyankar 
2775f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
2785f2c45f1SShri Abhyankar   header->ndata++;
2795f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2805f2c45f1SShri Abhyankar }
2815f2c45f1SShri Abhyankar 
2825f2c45f1SShri Abhyankar #undef __FUNCT__
2835f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetNumComponents"
2845f2c45f1SShri Abhyankar /*@
2855f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
2865f2c45f1SShri Abhyankar 
2875f2c45f1SShri Abhyankar   Not Collective
2885f2c45f1SShri Abhyankar 
2895f2c45f1SShri Abhyankar   Input Parameters:
2905f2c45f1SShri Abhyankar + dm - The DMNetwork object
2915f2c45f1SShri Abhyankar . p  - vertex/edge point
2925f2c45f1SShri Abhyankar 
2935f2c45f1SShri Abhyankar   Output Parameters:
2945f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
2955f2c45f1SShri Abhyankar 
2965f2c45f1SShri Abhyankar   Level: intermediate
2975f2c45f1SShri Abhyankar 
2985f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
2995f2c45f1SShri Abhyankar @*/
3005f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
3015f2c45f1SShri Abhyankar {
3025f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3035f2c45f1SShri Abhyankar   PetscInt       offset;
3045f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3055f2c45f1SShri Abhyankar 
3065f2c45f1SShri Abhyankar   PetscFunctionBegin;
3075f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
3085f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
3095f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3105f2c45f1SShri Abhyankar }
3115f2c45f1SShri Abhyankar 
3125f2c45f1SShri Abhyankar #undef __FUNCT__
3135f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetComponentTypeOffset"
3145f2c45f1SShri Abhyankar /*@
3155f2c45f1SShri Abhyankar   DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the
3165f2c45f1SShri Abhyankar                                     component value from the component data array
3175f2c45f1SShri Abhyankar 
3185f2c45f1SShri Abhyankar   Not Collective
3195f2c45f1SShri Abhyankar 
3205f2c45f1SShri Abhyankar   Input Parameters:
3215f2c45f1SShri Abhyankar + dm      - The DMNetwork object
3225f2c45f1SShri Abhyankar . p       - vertex/edge point
3235f2c45f1SShri Abhyankar - compnum - component number
3245f2c45f1SShri Abhyankar 
3255f2c45f1SShri Abhyankar   Output Parameters:
3265f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component
3275f2c45f1SShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
3285f2c45f1SShri Abhyankar 
3295f2c45f1SShri Abhyankar   Notes:
3305f2c45f1SShri Abhyankar   Typical usage:
3315f2c45f1SShri Abhyankar 
3325f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
3335f2c45f1SShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
3345f2c45f1SShri Abhyankar   Loop over vertices or edges
3355f2c45f1SShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
3365f2c45f1SShri Abhyankar     Loop over numcomps
3375f2c45f1SShri Abhyankar       DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset);
3385f2c45f1SShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
3395f2c45f1SShri Abhyankar 
3405f2c45f1SShri Abhyankar   Level: intermediate
3415f2c45f1SShri Abhyankar 
3425f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
3435f2c45f1SShri Abhyankar @*/
3445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
3455f2c45f1SShri Abhyankar {
3465f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
3475f2c45f1SShri Abhyankar   PetscInt                 offsetp;
3485f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
3495f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
3505f2c45f1SShri Abhyankar 
3515f2c45f1SShri Abhyankar   PetscFunctionBegin;
3525f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
3535f2c45f1SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
354d36f4e81SHong Zhang   if (compkey) *compkey = header->key[compnum];
355d36f4e81SHong Zhang   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
3565f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3575f2c45f1SShri Abhyankar }
3585f2c45f1SShri Abhyankar 
3595f2c45f1SShri Abhyankar #undef __FUNCT__
3605f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVariableOffset"
3615f2c45f1SShri Abhyankar /*@
3625f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
3635f2c45f1SShri Abhyankar 
3645f2c45f1SShri Abhyankar   Not Collective
3655f2c45f1SShri Abhyankar 
3665f2c45f1SShri Abhyankar   Input Parameters:
3675f2c45f1SShri Abhyankar + dm     - The DMNetwork object
3685f2c45f1SShri Abhyankar - p      - the edge/vertex point
3695f2c45f1SShri Abhyankar 
3705f2c45f1SShri Abhyankar   Output Parameters:
3715f2c45f1SShri Abhyankar . offset - the offset
3725f2c45f1SShri Abhyankar 
3735f2c45f1SShri Abhyankar   Level: intermediate
3745f2c45f1SShri Abhyankar 
3755f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
3765f2c45f1SShri Abhyankar @*/
3775f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
3785f2c45f1SShri Abhyankar {
3795f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3805f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3815f2c45f1SShri Abhyankar 
3825f2c45f1SShri Abhyankar   PetscFunctionBegin;
3835f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DofSection,p,offset);CHKERRQ(ierr);
3845f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3855f2c45f1SShri Abhyankar }
3865f2c45f1SShri Abhyankar 
3875f2c45f1SShri Abhyankar #undef __FUNCT__
3885f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVariableGlobalOffset"
3895f2c45f1SShri Abhyankar /*@
3905f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
3915f2c45f1SShri Abhyankar 
3925f2c45f1SShri Abhyankar   Not Collective
3935f2c45f1SShri Abhyankar 
3945f2c45f1SShri Abhyankar   Input Parameters:
3955f2c45f1SShri Abhyankar + dm      - The DMNetwork object
3965f2c45f1SShri Abhyankar - p       - the edge/vertex point
3975f2c45f1SShri Abhyankar 
3985f2c45f1SShri Abhyankar   Output Parameters:
3995f2c45f1SShri Abhyankar . offsetg - the offset
4005f2c45f1SShri Abhyankar 
4015f2c45f1SShri Abhyankar   Level: intermediate
4025f2c45f1SShri Abhyankar 
4035f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
4045f2c45f1SShri Abhyankar @*/
4055f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
4065f2c45f1SShri Abhyankar {
4075f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4085f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4095f2c45f1SShri Abhyankar 
4105f2c45f1SShri Abhyankar   PetscFunctionBegin;
4115f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);CHKERRQ(ierr);
412dd6f46cdSHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost node */
4135f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4145f2c45f1SShri Abhyankar }
4155f2c45f1SShri Abhyankar 
4165f2c45f1SShri Abhyankar #undef __FUNCT__
4175f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddNumVariables"
4185f2c45f1SShri Abhyankar /*@
4195f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
4205f2c45f1SShri Abhyankar 
4215f2c45f1SShri Abhyankar   Not Collective
4225f2c45f1SShri Abhyankar 
4235f2c45f1SShri Abhyankar   Input Parameters:
4245f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
4255f2c45f1SShri Abhyankar . p    - the vertex/edge point
4265f2c45f1SShri Abhyankar - nvar - number of additional variables
4275f2c45f1SShri Abhyankar 
4285f2c45f1SShri Abhyankar   Level: intermediate
4295f2c45f1SShri Abhyankar 
4305f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
4315f2c45f1SShri Abhyankar @*/
4325f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
4335f2c45f1SShri Abhyankar {
4345f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4355f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4365f2c45f1SShri Abhyankar 
4375f2c45f1SShri Abhyankar   PetscFunctionBegin;
4385f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
4395f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4405f2c45f1SShri Abhyankar }
4415f2c45f1SShri Abhyankar 
4425f2c45f1SShri Abhyankar #undef __FUNCT__
44327f51fceSHong Zhang #define __FUNCT__ "DMNetworkGetNumVariables"
44427f51fceSHong Zhang /*@
44527f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
44627f51fceSHong Zhang 
44727f51fceSHong Zhang   Not Collective
44827f51fceSHong Zhang 
44927f51fceSHong Zhang   Input Parameters:
45027f51fceSHong Zhang + dm   - The DMNetworkObject
45127f51fceSHong Zhang - p    - the vertex/edge point
45227f51fceSHong Zhang 
45327f51fceSHong Zhang   Output Parameters:
45427f51fceSHong Zhang . nvar - number of variables
45527f51fceSHong Zhang 
45627f51fceSHong Zhang   Level: intermediate
45727f51fceSHong Zhang 
45827f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
45927f51fceSHong Zhang @*/
46027f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
46127f51fceSHong Zhang {
46227f51fceSHong Zhang   PetscErrorCode ierr;
46327f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
46427f51fceSHong Zhang 
46527f51fceSHong Zhang   PetscFunctionBegin;
46627f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
46727f51fceSHong Zhang   PetscFunctionReturn(0);
46827f51fceSHong Zhang }
46927f51fceSHong Zhang 
47027f51fceSHong Zhang #undef __FUNCT__
471662b5b01SHong Zhang #define __FUNCT__ "DMNetworkSetNumVariables"
4725f2c45f1SShri Abhyankar /*@
4735f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
4745f2c45f1SShri Abhyankar 
4755f2c45f1SShri Abhyankar   Not Collective
4765f2c45f1SShri Abhyankar 
4775f2c45f1SShri Abhyankar   Input Parameters:
4785f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
4795f2c45f1SShri Abhyankar . p    - the vertex/edge point
4805f2c45f1SShri Abhyankar - nvar - number of variables
4815f2c45f1SShri Abhyankar 
4825f2c45f1SShri Abhyankar   Level: intermediate
4835f2c45f1SShri Abhyankar 
4845f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
4855f2c45f1SShri Abhyankar @*/
4865f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
4875f2c45f1SShri Abhyankar {
4885f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4895f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4905f2c45f1SShri Abhyankar 
4915f2c45f1SShri Abhyankar   PetscFunctionBegin;
4925f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
4935f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4945f2c45f1SShri Abhyankar }
4955f2c45f1SShri Abhyankar 
4965f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
4975f2c45f1SShri Abhyankar    function is called during DMSetUp() */
4985f2c45f1SShri Abhyankar #undef __FUNCT__
4995f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkComponentSetUp"
5005f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
5015f2c45f1SShri Abhyankar {
5025f2c45f1SShri Abhyankar   PetscErrorCode              ierr;
5035f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5045f2c45f1SShri Abhyankar   PetscInt                    arr_size;
5055f2c45f1SShri Abhyankar   PetscInt                    p,offset,offsetp;
5065f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
5075f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
5085f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType      *componentdataarray;
5095f2c45f1SShri Abhyankar   PetscInt ncomp, i;
5105f2c45f1SShri Abhyankar 
5115f2c45f1SShri Abhyankar   PetscFunctionBegin;
5125f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
5135f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
51475b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
5155f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
5165f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
5175f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
5185f2c45f1SShri Abhyankar     /* Copy header */
5195f2c45f1SShri Abhyankar     header = &network->header[p];
520302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
5215f2c45f1SShri Abhyankar     /* Copy data */
5225f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
5235f2c45f1SShri Abhyankar     ncomp = header->ndata;
5245f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
5255f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
526302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
5275f2c45f1SShri Abhyankar     }
5285f2c45f1SShri Abhyankar   }
5295f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5305f2c45f1SShri Abhyankar }
5315f2c45f1SShri Abhyankar 
5325f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
5335f2c45f1SShri Abhyankar #undef __FUNCT__
5345f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkVariablesSetUp"
5355f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
5365f2c45f1SShri Abhyankar {
5375f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5385f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5395f2c45f1SShri Abhyankar 
5405f2c45f1SShri Abhyankar   PetscFunctionBegin;
5415f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
5425f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5435f2c45f1SShri Abhyankar }
5445f2c45f1SShri Abhyankar 
5455f2c45f1SShri Abhyankar #undef __FUNCT__
5465f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetComponentDataArray"
5475f2c45f1SShri Abhyankar /*@C
5485f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
5495f2c45f1SShri Abhyankar 
5505f2c45f1SShri Abhyankar   Not Collective
5515f2c45f1SShri Abhyankar 
5525f2c45f1SShri Abhyankar   Input Parameters:
5535f2c45f1SShri Abhyankar . dm - The DMNetwork Object
5545f2c45f1SShri Abhyankar 
5555f2c45f1SShri Abhyankar   Output Parameters:
5565f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
5575f2c45f1SShri Abhyankar 
5585f2c45f1SShri Abhyankar   Level: intermediate
5595f2c45f1SShri Abhyankar 
5605f2c45f1SShri Abhyankar .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents
5615f2c45f1SShri Abhyankar @*/
5625f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
5635f2c45f1SShri Abhyankar {
5645f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5655f2c45f1SShri Abhyankar 
5665f2c45f1SShri Abhyankar   PetscFunctionBegin;
5675f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
5685f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5695f2c45f1SShri Abhyankar }
5705f2c45f1SShri Abhyankar 
5715f2c45f1SShri Abhyankar #undef __FUNCT__
5725f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkDistribute"
5735f2c45f1SShri Abhyankar /*@
5745f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
5755f2c45f1SShri Abhyankar 
5765f2c45f1SShri Abhyankar   Collective
5775f2c45f1SShri Abhyankar 
5785f2c45f1SShri Abhyankar   Input Parameter:
5795f2c45f1SShri Abhyankar + oldDM - the original DMNetwork object
5805f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
5815f2c45f1SShri Abhyankar 
5825f2c45f1SShri Abhyankar   Output Parameter:
5835f2c45f1SShri Abhyankar . distDM - the distributed DMNetwork object
5845f2c45f1SShri Abhyankar 
5855f2c45f1SShri Abhyankar   Notes:
5865f2c45f1SShri Abhyankar   This routine should be called only when using multiple processors.
5875f2c45f1SShri Abhyankar 
5888b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
5895f2c45f1SShri Abhyankar 
5905f2c45f1SShri Abhyankar   Level: intermediate
5915f2c45f1SShri Abhyankar 
5925f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
5935f2c45f1SShri Abhyankar @*/
59480cf41d5SMatthew G. Knepley PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM)
5955f2c45f1SShri Abhyankar {
5965f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5975f2c45f1SShri Abhyankar   DM_Network     *oldDMnetwork = (DM_Network*)oldDM->data;
5985f2c45f1SShri Abhyankar   PetscSF        pointsf;
5995f2c45f1SShri Abhyankar   DM             newDM;
6005f2c45f1SShri Abhyankar   DM_Network     *newDMnetwork;
6015f2c45f1SShri Abhyankar 
6025f2c45f1SShri Abhyankar   PetscFunctionBegin;
6035f2c45f1SShri Abhyankar   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr);
6045f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
6055f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
6065f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
60780cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
6085f2c45f1SShri Abhyankar   /* Distribute dof section */
6095f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr);
6105f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
6115f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr);
6125f2c45f1SShri Abhyankar   /* Distribute data and associated section */
61331da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
6145f2c45f1SShri Abhyankar   /* Destroy point SF */
6155f2c45f1SShri Abhyankar   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
6165f2c45f1SShri Abhyankar 
6175f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
6185f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
6195f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
6205f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
6215f2c45f1SShri Abhyankar   newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart;
6225f2c45f1SShri Abhyankar   newDMnetwork->NNodes = oldDMnetwork->NNodes;
6235f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
6245f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
6255f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
6265f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
6275f2c45f1SShri Abhyankar 
6285f2c45f1SShri Abhyankar   *distDM = newDM;
6295f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6305f2c45f1SShri Abhyankar }
6315f2c45f1SShri Abhyankar 
6325f2c45f1SShri Abhyankar #undef __FUNCT__
6335f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetSupportingEdges"
6345f2c45f1SShri Abhyankar /*@C
6355f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
6365f2c45f1SShri Abhyankar 
6375f2c45f1SShri Abhyankar   Not Collective
6385f2c45f1SShri Abhyankar 
6395f2c45f1SShri Abhyankar   Input Parameters:
6405f2c45f1SShri Abhyankar + dm - The DMNetwork object
6415f2c45f1SShri Abhyankar - p  - the vertex point
6425f2c45f1SShri Abhyankar 
6435f2c45f1SShri Abhyankar   Output Paramters:
6445f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
6455f2c45f1SShri Abhyankar - edges  - List of edge points
6465f2c45f1SShri Abhyankar 
6475f2c45f1SShri Abhyankar   Level: intermediate
6485f2c45f1SShri Abhyankar 
6495f2c45f1SShri Abhyankar   Fortran Notes:
6505f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
6515f2c45f1SShri Abhyankar   include petsc.h90 in your code.
6525f2c45f1SShri Abhyankar 
6535f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
6545f2c45f1SShri Abhyankar @*/
6555f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
6565f2c45f1SShri Abhyankar {
6575f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6585f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6595f2c45f1SShri Abhyankar 
6605f2c45f1SShri Abhyankar   PetscFunctionBegin;
6615f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
6625f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
6635f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6645f2c45f1SShri Abhyankar }
6655f2c45f1SShri Abhyankar 
6665f2c45f1SShri Abhyankar #undef __FUNCT__
6675f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetConnectedNodes"
6685f2c45f1SShri Abhyankar /*@C
669596e729fSHong Zhang   DMNetworkGetConnectedNodes - Return the connected vertices for this edge point
6705f2c45f1SShri Abhyankar 
6715f2c45f1SShri Abhyankar   Not Collective
6725f2c45f1SShri Abhyankar 
6735f2c45f1SShri Abhyankar   Input Parameters:
6745f2c45f1SShri Abhyankar + dm - The DMNetwork object
6755f2c45f1SShri Abhyankar - p  - the edge point
6765f2c45f1SShri Abhyankar 
6775f2c45f1SShri Abhyankar   Output Paramters:
6785f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
6795f2c45f1SShri Abhyankar 
6805f2c45f1SShri Abhyankar   Level: intermediate
6815f2c45f1SShri Abhyankar 
6825f2c45f1SShri Abhyankar   Fortran Notes:
6835f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
6845f2c45f1SShri Abhyankar   include petsc.h90 in your code.
6855f2c45f1SShri Abhyankar 
6865f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
6875f2c45f1SShri Abhyankar @*/
6885f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
6895f2c45f1SShri Abhyankar {
6905f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6915f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6925f2c45f1SShri Abhyankar 
6935f2c45f1SShri Abhyankar   PetscFunctionBegin;
6945f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
6955f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6965f2c45f1SShri Abhyankar }
6975f2c45f1SShri Abhyankar 
6985f2c45f1SShri Abhyankar #undef __FUNCT__
6995f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkIsGhostVertex"
7005f2c45f1SShri Abhyankar /*@
7015f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
7025f2c45f1SShri Abhyankar 
7035f2c45f1SShri Abhyankar   Not Collective
7045f2c45f1SShri Abhyankar 
7055f2c45f1SShri Abhyankar   Input Parameters:
7065f2c45f1SShri Abhyankar + dm - The DMNetwork object
7075f2c45f1SShri Abhyankar . p  - the vertex point
7085f2c45f1SShri Abhyankar 
7095f2c45f1SShri Abhyankar   Output Parameter:
7105f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
7115f2c45f1SShri Abhyankar 
7125f2c45f1SShri Abhyankar   Level: intermediate
7135f2c45f1SShri Abhyankar 
7145f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
7155f2c45f1SShri Abhyankar @*/
7165f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
7175f2c45f1SShri Abhyankar {
7185f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7195f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
7205f2c45f1SShri Abhyankar   PetscInt       offsetg;
7215f2c45f1SShri Abhyankar   PetscSection   sectiong;
7225f2c45f1SShri Abhyankar 
7235f2c45f1SShri Abhyankar   PetscFunctionBegin;
7245f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
7255f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
7265f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
7275f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
7285f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7295f2c45f1SShri Abhyankar }
7305f2c45f1SShri Abhyankar 
7315f2c45f1SShri Abhyankar #undef __FUNCT__
7325f2c45f1SShri Abhyankar #define __FUNCT__ "DMSetUp_Network"
7335f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
7345f2c45f1SShri Abhyankar {
7355f2c45f1SShri Abhyankar   PetscErrorCode ierr;
7365f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
7375f2c45f1SShri Abhyankar 
7385f2c45f1SShri Abhyankar   PetscFunctionBegin;
7395f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
7405f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
7415f2c45f1SShri Abhyankar 
7425f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
7435f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
7445f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
7455f2c45f1SShri Abhyankar }
7465f2c45f1SShri Abhyankar 
7475f2c45f1SShri Abhyankar #undef __FUNCT__
74817df6e9eSHong Zhang #define __FUNCT__ "DMNetworkHasJacobian"
7491ad426b7SHong Zhang /*@
75017df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
7511ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
7521ad426b7SHong Zhang 
7531ad426b7SHong Zhang     Collective
7541ad426b7SHong Zhang 
7551ad426b7SHong Zhang     Input Parameters:
75683b2e829SHong Zhang +   dm - The DMNetwork object
75783b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
75883b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
7591ad426b7SHong Zhang 
7601ad426b7SHong Zhang     Level: intermediate
7611ad426b7SHong Zhang 
7621ad426b7SHong Zhang @*/
76383b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
7641ad426b7SHong Zhang {
7651ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
7661ad426b7SHong Zhang 
7671ad426b7SHong Zhang   PetscFunctionBegin;
76883b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
76983b2e829SHong Zhang   network->userVertexJacobian = vflg;
7701ad426b7SHong Zhang   PetscFunctionReturn(0);
7711ad426b7SHong Zhang }
7721ad426b7SHong Zhang 
7731ad426b7SHong Zhang #undef __FUNCT__
77483b2e829SHong Zhang #define __FUNCT__ "DMNetworkEdgeSetMatrix"
7751ad426b7SHong Zhang /*@
77683b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
77783b2e829SHong Zhang 
77883b2e829SHong Zhang     Not Collective
77983b2e829SHong Zhang 
78083b2e829SHong Zhang     Input Parameters:
78183b2e829SHong Zhang +   dm - The DMNetwork object
78283b2e829SHong Zhang .   p  - the edge point
7833e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
7843e97b6e8SHong Zhang         J[0]: this edge
7853e97b6e8SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes()
78683b2e829SHong Zhang 
78783b2e829SHong Zhang     Level: intermediate
78883b2e829SHong Zhang 
78983b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
79083b2e829SHong Zhang @*/
79183b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
79283b2e829SHong Zhang {
79383b2e829SHong Zhang   PetscErrorCode ierr;
79483b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
79583b2e829SHong Zhang 
79683b2e829SHong Zhang   PetscFunctionBegin;
797883e35e8SHong Zhang   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
79883b2e829SHong Zhang   if (!network->Je) {
799883e35e8SHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
80083b2e829SHong Zhang   }
801883e35e8SHong Zhang   network->Je[3*p]   = J[0];
802883e35e8SHong Zhang   network->Je[3*p+1] = J[1];
803883e35e8SHong Zhang   network->Je[3*p+2] = J[2];
80483b2e829SHong Zhang   PetscFunctionReturn(0);
80583b2e829SHong Zhang }
80683b2e829SHong Zhang 
80783b2e829SHong Zhang #undef __FUNCT__
808883e35e8SHong Zhang #define __FUNCT__ "DMNetworkVertexSetMatrix"
80983b2e829SHong Zhang /*@
81076ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
8111ad426b7SHong Zhang 
8121ad426b7SHong Zhang     Not Collective
8131ad426b7SHong Zhang 
8141ad426b7SHong Zhang     Input Parameters:
8151ad426b7SHong Zhang +   dm - The DMNetwork object
8161ad426b7SHong Zhang .   p  - the vertex point
8173e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
8183e97b6e8SHong Zhang         J[0]:       this vertex
8193e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
8203e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
8211ad426b7SHong Zhang 
8221ad426b7SHong Zhang     Level: intermediate
8231ad426b7SHong Zhang 
82483b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
8251ad426b7SHong Zhang @*/
826883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
8275f2c45f1SShri Abhyankar {
8285f2c45f1SShri Abhyankar   PetscErrorCode ierr;
8295f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
830f4431b8cSHong Zhang   PetscInt       i,*vptr,nedges,vStart,vEnd;
831883e35e8SHong Zhang   const PetscInt *edges;
8325f2c45f1SShri Abhyankar 
8335f2c45f1SShri Abhyankar   PetscFunctionBegin;
834883e35e8SHong Zhang   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
835883e35e8SHong Zhang 
836883e35e8SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
837f4431b8cSHong Zhang 
83883b2e829SHong Zhang   if (!network->Jv) {
839f4431b8cSHong Zhang     PetscInt nNodes = network->nNodes,nedges_total;
840883e35e8SHong Zhang     /* count nvertex_total */
8413e97b6e8SHong Zhang     nedges_total = 0;
842883e35e8SHong Zhang     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
843f4431b8cSHong Zhang     ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr);
844f4431b8cSHong Zhang 
845883e35e8SHong Zhang     vptr[0] = 0;
846f4431b8cSHong Zhang     for (i=0; i<nNodes; i++) {
847f4431b8cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
848883e35e8SHong Zhang       nedges_total += nedges;
849f4431b8cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
8501ad426b7SHong Zhang     }
8513e97b6e8SHong Zhang 
852f4431b8cSHong Zhang     ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr);
853883e35e8SHong Zhang     network->Jvptr = vptr;
854883e35e8SHong Zhang   }
855883e35e8SHong Zhang 
856883e35e8SHong Zhang   vptr = network->Jvptr;
8573e97b6e8SHong Zhang   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
8583e97b6e8SHong Zhang 
8593e97b6e8SHong Zhang   /* Set Jacobian for each supporting edge and connected vertex */
860883e35e8SHong Zhang   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
861883e35e8SHong Zhang   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
862883e35e8SHong Zhang   PetscFunctionReturn(0);
863883e35e8SHong Zhang }
864883e35e8SHong Zhang 
865883e35e8SHong Zhang #undef __FUNCT__
866883e35e8SHong Zhang #define __FUNCT__ "MatSetDenseblock_private"
8673e97b6e8SHong Zhang PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
868883e35e8SHong Zhang {
869883e35e8SHong Zhang   PetscErrorCode ierr;
870883e35e8SHong Zhang   PetscInt       j,*cols;
871883e35e8SHong Zhang   PetscScalar    *zeros;
872883e35e8SHong Zhang 
873883e35e8SHong Zhang   PetscFunctionBegin;
874883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
875883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
876883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
877883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
8781ad426b7SHong Zhang   PetscFunctionReturn(0);
8791ad426b7SHong Zhang }
880a4e85ca8SHong Zhang 
8813e97b6e8SHong Zhang #undef __FUNCT__
8823e97b6e8SHong Zhang #define __FUNCT__ "MatSetUserblock_private"
8833e97b6e8SHong Zhang PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
8843e97b6e8SHong Zhang {
8853e97b6e8SHong Zhang   PetscErrorCode ierr;
8863e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
8873e97b6e8SHong Zhang   const PetscInt *cols;
8883e97b6e8SHong Zhang   PetscScalar    zero=0.0;
8893e97b6e8SHong Zhang 
8903e97b6e8SHong Zhang   PetscFunctionBegin;
8913e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
8923e97b6e8SHong Zhang   if (nrows != M || ncols != N) SETERRQ4(PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N);
8933e97b6e8SHong Zhang 
8943e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
8953e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
8963e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
8973e97b6e8SHong Zhang       col = cols[j] + cstart;
8983e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
8993e97b6e8SHong Zhang     }
9003e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
9013e97b6e8SHong Zhang   }
9023e97b6e8SHong Zhang   PetscFunctionReturn(0);
9033e97b6e8SHong Zhang }
9041ad426b7SHong Zhang 
9051ad426b7SHong Zhang #undef __FUNCT__
906a4e85ca8SHong Zhang #define __FUNCT__ "MatSetblock_private"
907a4e85ca8SHong Zhang PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
908a4e85ca8SHong Zhang {
909a4e85ca8SHong Zhang   PetscErrorCode ierr;
910f4431b8cSHong Zhang 
911a4e85ca8SHong Zhang   PetscFunctionBegin;
912a4e85ca8SHong Zhang   if (Ju) {
913a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
914a4e85ca8SHong Zhang   } else {
915a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
916a4e85ca8SHong Zhang   }
917a4e85ca8SHong Zhang   PetscFunctionReturn(0);
918a4e85ca8SHong Zhang }
919a4e85ca8SHong Zhang 
920a4e85ca8SHong Zhang #undef __FUNCT__
9211ad426b7SHong Zhang #define __FUNCT__ "DMCreateMatrix_Network"
9221ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
9231ad426b7SHong Zhang {
9241ad426b7SHong Zhang   PetscErrorCode ierr;
9251ad426b7SHong Zhang   DM_Network     *network = (DM_Network*) dm->data;
926a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
927*840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
9281ad426b7SHong Zhang   PetscBool      ghost;
929a4e85ca8SHong Zhang   Mat            Juser;
930bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
931a4e85ca8SHong Zhang   PetscInt       nedges,*vptr,vc;
932a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
9331ad426b7SHong Zhang 
9341ad426b7SHong Zhang   PetscFunctionBegin;
935a4e85ca8SHong Zhang   if (!network->userEdgeJacobian && !network->userVertexJacobian) {
936a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
937bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
938bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
9391ad426b7SHong Zhang     PetscFunctionReturn(0);
9401ad426b7SHong Zhang   }
9411ad426b7SHong Zhang 
942bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
9432a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
944bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
945bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
9462a945128SHong Zhang 
9472a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
9482a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
949*840c2264SHong Zhang #if 0
950a4e85ca8SHong Zhang   /* Preallocation - submatrix for an element (edge/vertex) is allocated as a dense block,
951a4e85ca8SHong Zhang    see DMCreateMatrix_Plex() */
952*840c2264SHong Zhang   PetscInt *dnz,*onz,*dnzu,*onzu;
953bfbc38dcSHong Zhang   ierr = PetscCalloc4(localSize,&dnz,localSize,&onz,localSize,&dnzu,localSize,&onzu);CHKERRQ(ierr);
954bfbc38dcSHong Zhang   ierr = DMPlexPreallocateOperator(network->plex,1,dnz,onz,dnzu,onzu,*J,PETSC_FALSE);CHKERRQ(ierr);
955bfbc38dcSHong Zhang   ierr = PetscFree4(dnz,onz,dnzu,onzu);CHKERRQ(ierr);
956*840c2264SHong Zhang #endif
957*840c2264SHong Zhang   //========================= Sparse preallocation =======================================
958*840c2264SHong Zhang   MPI_Comm    comm;
959*840c2264SHong Zhang   PetscMPIInt rank;
960*840c2264SHong Zhang   Vec         vd_nz,vo_nz;
961*840c2264SHong Zhang   PetscInt    M;
962*840c2264SHong Zhang   PetscScalar val;
963*840c2264SHong Zhang 
964*840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
965*840c2264SHong Zhang   ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
966*840c2264SHong Zhang 
967*840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
968*840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
969*840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
970*840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
971*840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
972*840c2264SHong Zhang 
973*840c2264SHong Zhang   ierr = VecGetSize(vd_nz,&M);CHKERRQ(ierr);
974*840c2264SHong Zhang   printf("[%d] J localSize %d, gSize %d\n",rank,localSize,M);
975*840c2264SHong Zhang 
976*840c2264SHong Zhang   /* Set matrix entries for edges */
977*840c2264SHong Zhang   /*------------------------------*/
978*840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
979*840c2264SHong Zhang 
980*840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
981*840c2264SHong Zhang     /* Get row indices */
982*840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
983*840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
984*840c2264SHong Zhang     if (nrows) {
985*840c2264SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
986*840c2264SHong Zhang       ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
987*840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
988*840c2264SHong Zhang 
989*840c2264SHong Zhang       /* Set matrix entries for conntected vertices */
990*840c2264SHong Zhang       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
991*840c2264SHong Zhang       for (v=0; v<2; v++) {
992*840c2264SHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
993*840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
994*840c2264SHong Zhang 
995*840c2264SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
996*840c2264SHong Zhang         //ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
997*840c2264SHong Zhang         val = (PetscScalar)ncols;
998*840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
999*840c2264SHong Zhang         if (!ghost) {
1000*840c2264SHong Zhang           //printf("[%d] e%d - v%d ncols %d\n",rank,e,cone[v],ncols);
1001*840c2264SHong Zhang           for (j=0; j<nrows; j++) {
1002*840c2264SHong Zhang             ierr = VecSetValues(vd_nz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1003*840c2264SHong Zhang           }
1004*840c2264SHong Zhang         } else {
1005*840c2264SHong Zhang           //printf("[%d] e%d - v%d is ghost, ncols %d\n",rank,e,cone[v],ncols);
1006*840c2264SHong Zhang           for (j=0; j<nrows; j++) {
1007*840c2264SHong Zhang             ierr = VecSetValues(vo_nz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1008*840c2264SHong Zhang           }
1009*840c2264SHong Zhang         }
1010*840c2264SHong Zhang       }
1011*840c2264SHong Zhang 
1012*840c2264SHong Zhang       /* Set matrix entries for edge self */
1013*840c2264SHong Zhang       cstart = rstart;
1014*840c2264SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1015*840c2264SHong Zhang 
1016*840c2264SHong Zhang       //ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1017*840c2264SHong Zhang       //printf("[%d] e%d, vd_nz %d\n",rank,e,nrows);
1018*840c2264SHong Zhang       for (j=0; j<nrows; j++) {
1019*840c2264SHong Zhang         PetscInt ncols_u;
1020*840c2264SHong Zhang         ierr = MatGetRow(Juser,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1021*840c2264SHong Zhang         val = (PetscScalar)ncols_u;
1022*840c2264SHong Zhang         ierr = VecSetValues(vd_nz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1023*840c2264SHong Zhang         ierr = MatRestoreRow(Juser,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
1024*840c2264SHong Zhang       }
1025*840c2264SHong Zhang       ierr = PetscFree(rows);CHKERRQ(ierr);
1026*840c2264SHong Zhang     }
1027*840c2264SHong Zhang   }
1028*840c2264SHong Zhang 
1029*840c2264SHong Zhang   /* Set matrix entries for vertices */
1030*840c2264SHong Zhang   /*---------------------------------*/
1031*840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1032*840c2264SHong Zhang   if (vEnd - vStart) {
1033*840c2264SHong Zhang     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1034*840c2264SHong Zhang     vptr = network->Jvptr;
1035*840c2264SHong Zhang     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1036*840c2264SHong Zhang   }
1037*840c2264SHong Zhang 
1038*840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1039*840c2264SHong Zhang     /* Get row indices */
1040*840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1041*840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1042*840c2264SHong Zhang     if (!nrows) continue;
1043*840c2264SHong Zhang 
1044*840c2264SHong Zhang     ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
1045*840c2264SHong Zhang     for (j=0; j<nrows; j++) rows[j] = j + rstart;
1046*840c2264SHong Zhang 
1047*840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1048*840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1049*840c2264SHong Zhang 
1050*840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1051*840c2264SHong Zhang     //if (ghost) printf("[%d] v%d is ghost, nedges %d\n",rank,v,nedges);
1052*840c2264SHong Zhang 
1053*840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1054*840c2264SHong Zhang       /* Supporting edges */
1055*840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1056*840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1057*840c2264SHong Zhang 
1058*840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1059*840c2264SHong Zhang       //ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1060*840c2264SHong Zhang       val = (PetscScalar)ncols;
1061*840c2264SHong Zhang       if (!ghost) {
1062*840c2264SHong Zhang         //printf("[%d] v%d supporting e%d ncols %d\n",rank,v,edges[e],ncols);
1063*840c2264SHong Zhang         for (j=0; j<nrows; j++) {
1064*840c2264SHong Zhang           ierr = VecSetValues(vd_nz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1065*840c2264SHong Zhang         }
1066*840c2264SHong Zhang       } else {
1067*840c2264SHong Zhang         //printf("[%d] v%d is ghost, supporting e%d ncols %d\n",rank,v,edges[e],ncols);
1068*840c2264SHong Zhang         for (j=0; j<nrows; j++) {
1069*840c2264SHong Zhang           ierr = VecSetValues(vo_nz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1070*840c2264SHong Zhang         }
1071*840c2264SHong Zhang       }
1072*840c2264SHong Zhang 
1073*840c2264SHong Zhang       /* Connected vertices */
1074*840c2264SHong Zhang       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1075*840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1076*840c2264SHong Zhang       PetscBool ghost_vc;
1077*840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1078*840c2264SHong Zhang       //if (ghost_vc) printf("[%d] gost v%d is connected to v%d\n",rank,v,vc);
1079*840c2264SHong Zhang 
1080*840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
1081*840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1082*840c2264SHong Zhang 
1083*840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1084*840c2264SHong Zhang       //ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1085*840c2264SHong Zhang       val = (PetscScalar)ncols;
1086*840c2264SHong Zhang       if (!ghost_vc && !ghost) {
1087*840c2264SHong Zhang         //printf("[%d] v%d connected vc%d ncols %d\n",rank,v,vc,ncols);
1088*840c2264SHong Zhang         for (j=0; j<nrows; j++) {
1089*840c2264SHong Zhang           //vals[j] = (PetscScalar)ncols;
1090*840c2264SHong Zhang           ierr = VecSetValues(vd_nz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1091*840c2264SHong Zhang         }
1092*840c2264SHong Zhang       } else {
1093*840c2264SHong Zhang         //printf("[%d] v%d is ghost, connected vc%d ncols %d\n",rank,v,vc,ncols);
1094*840c2264SHong Zhang         for (j=0; j<nrows; j++) {
1095*840c2264SHong Zhang           //vals[j] = (PetscScalar)ncols;
1096*840c2264SHong Zhang           ierr = VecSetValues(vo_nz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1097*840c2264SHong Zhang         }
1098*840c2264SHong Zhang       }
1099*840c2264SHong Zhang     }
1100*840c2264SHong Zhang 
1101*840c2264SHong Zhang     /* Set matrix entries for vertex self */
1102*840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1103*840c2264SHong Zhang     if (!ghost) {
1104*840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1105*840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1106*840c2264SHong Zhang 
1107*840c2264SHong Zhang       //ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1108*840c2264SHong Zhang       //printf("[%d] v%d ncols %d\n",rank,v,ncols);
1109*840c2264SHong Zhang       val = (PetscScalar)nrows;
1110*840c2264SHong Zhang       for (j=0; j<nrows; j++) {
1111*840c2264SHong Zhang         ierr = VecSetValues(vd_nz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
1112*840c2264SHong Zhang       }
1113*840c2264SHong Zhang     }
1114*840c2264SHong Zhang     ierr = PetscFree(rows);CHKERRQ(ierr);
1115*840c2264SHong Zhang   }
1116*840c2264SHong Zhang 
1117*840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1118*840c2264SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1119*840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
1120*840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1121*840c2264SHong Zhang 
1122*840c2264SHong Zhang   PetscInt    dnnz[localSize],onnz[localSize];
1123*840c2264SHong Zhang   PetscScalar *vdnz,*vonz;
1124*840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1125*840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1126*840c2264SHong Zhang 
1127*840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1128*840c2264SHong Zhang     dnnz[j] = (PetscInt)vdnz[j];
1129*840c2264SHong Zhang     onnz[j] = (PetscInt)vonz[j];
1130*840c2264SHong Zhang   }
1131*840c2264SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
1132*840c2264SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
1133*840c2264SHong Zhang #if 0
1134*840c2264SHong Zhang   if (!rank) printf("vd_nz\n");
1135*840c2264SHong Zhang   ierr = VecView(vd_nz,0);CHKERRQ(ierr);
1136*840c2264SHong Zhang   if (!rank) printf("vo_nz\n");
1137*840c2264SHong Zhang   ierr = VecView(vo_nz,0);CHKERRQ(ierr);
1138*840c2264SHong Zhang #endif
1139*840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1140*840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1141*840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1142*840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1143*840c2264SHong Zhang 
1144*840c2264SHong Zhang   //=====================================  endof new
1145bfbc38dcSHong Zhang 
1146bfbc38dcSHong Zhang   /* Set matrix entries for edges */
114783b2e829SHong Zhang   /*------------------------------*/
1148bfbc38dcSHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
11494b976069SHong Zhang 
11501ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1151bfbc38dcSHong Zhang     /* Get row indices */
11521ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
115317df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
11544b976069SHong Zhang     if (nrows) {
11554b976069SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1156a4e85ca8SHong Zhang       ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
115717df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
11581ad426b7SHong Zhang 
1159bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
116017df6e9eSHong Zhang       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1161bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1162bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1163883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
11643e97b6e8SHong Zhang 
1165a4e85ca8SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1166a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1167bfbc38dcSHong Zhang       }
116817df6e9eSHong Zhang 
1169bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
11703e97b6e8SHong Zhang       cstart = rstart;
1171a4e85ca8SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1172a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
1173a4e85ca8SHong Zhang       ierr = PetscFree(rows);CHKERRQ(ierr);
11741ad426b7SHong Zhang     }
11754b976069SHong Zhang   }
11761ad426b7SHong Zhang 
1177bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
117883b2e829SHong Zhang   /*---------------------------------*/
1179bfbc38dcSHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
11804b976069SHong Zhang   if (vEnd - vStart) {
11812a808120SBarry Smith     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"User must provide Jv");
11822a808120SBarry Smith     if (!(vptr = network->Jvptr)) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"User must provide vptr");
11834b976069SHong Zhang 
11841ad426b7SHong Zhang     for (v=vStart; v<vEnd; v++) {
1185bfbc38dcSHong Zhang       /* Get row indices */
1186596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1187596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
11884b976069SHong Zhang       if (!nrows) continue;
1189596e729fSHong Zhang 
1190a4e85ca8SHong Zhang       ierr = PetscMalloc1(nrows,&rows);CHKERRQ(ierr);
1191596e729fSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1192596e729fSHong Zhang 
1193bfbc38dcSHong Zhang       /* Get supporting edges and connected vertices */
1194596e729fSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1195596e729fSHong Zhang 
1196596e729fSHong Zhang       for (e=0; e<nedges; e++) {
1197bfbc38dcSHong Zhang         /* Supporting edges */
1198596e729fSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1199596e729fSHong Zhang         ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1200596e729fSHong Zhang 
1201a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1202a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1203596e729fSHong Zhang 
1204bfbc38dcSHong Zhang         /* Connected vertices */
120544aca652SHong Zhang         ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
12062a945128SHong Zhang         vc = (v == cone[0]) ? cone[1]:cone[0];
12072a945128SHong Zhang 
120844aca652SHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
120944aca652SHong Zhang         ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1210a4e85ca8SHong Zhang 
1211a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1212*840c2264SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1213596e729fSHong Zhang       }
1214596e729fSHong Zhang 
1215bfbc38dcSHong Zhang       /* Set matrix entries for vertex self */
12161ad426b7SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
12171ad426b7SHong Zhang       if (!ghost) {
1218596e729fSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1219a4e85ca8SHong Zhang         Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1220a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
12211ad426b7SHong Zhang       }
1222a4e85ca8SHong Zhang       ierr = PetscFree(rows);CHKERRQ(ierr);
12231ad426b7SHong Zhang     }
1224c619b0b1SBarry Smith   }
12251ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
12261ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1227dd6f46cdSHong Zhang 
12285f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
12295f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12305f2c45f1SShri Abhyankar }
12315f2c45f1SShri Abhyankar 
12325f2c45f1SShri Abhyankar #undef __FUNCT__
12335f2c45f1SShri Abhyankar #define __FUNCT__ "DMDestroy_Network"
12345f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
12355f2c45f1SShri Abhyankar {
12365f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12375f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
12385f2c45f1SShri Abhyankar 
12395f2c45f1SShri Abhyankar   PetscFunctionBegin;
12408415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
124183b2e829SHong Zhang   if (network->Je) {
124283b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
124383b2e829SHong Zhang   }
124483b2e829SHong Zhang   if (network->Jv) {
1245883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
124683b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
12471ad426b7SHong Zhang   }
12485f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
12495f2c45f1SShri Abhyankar   network->edges = NULL;
12505f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
12515f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
125283b2e829SHong Zhang 
12535f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
12545f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
12555f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
12565f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
12575f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12585f2c45f1SShri Abhyankar }
12595f2c45f1SShri Abhyankar 
12605f2c45f1SShri Abhyankar #undef __FUNCT__
12615f2c45f1SShri Abhyankar #define __FUNCT__ "DMView_Network"
12625f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
12635f2c45f1SShri Abhyankar {
12645f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12655f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
12665f2c45f1SShri Abhyankar 
12675f2c45f1SShri Abhyankar   PetscFunctionBegin;
12685f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
12695f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12705f2c45f1SShri Abhyankar }
12715f2c45f1SShri Abhyankar 
12725f2c45f1SShri Abhyankar #undef __FUNCT__
12735f2c45f1SShri Abhyankar #define __FUNCT__ "DMGlobalToLocalBegin_Network"
12745f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
12755f2c45f1SShri Abhyankar {
12765f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12775f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
12785f2c45f1SShri Abhyankar 
12795f2c45f1SShri Abhyankar   PetscFunctionBegin;
12805f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
12815f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12825f2c45f1SShri Abhyankar }
12835f2c45f1SShri Abhyankar 
12845f2c45f1SShri Abhyankar #undef __FUNCT__
12855f2c45f1SShri Abhyankar #define __FUNCT__ "DMGlobalToLocalEnd_Network"
12865f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
12875f2c45f1SShri Abhyankar {
12885f2c45f1SShri Abhyankar   PetscErrorCode ierr;
12895f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
12905f2c45f1SShri Abhyankar 
12915f2c45f1SShri Abhyankar   PetscFunctionBegin;
12925f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
12935f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
12945f2c45f1SShri Abhyankar }
12955f2c45f1SShri Abhyankar 
12965f2c45f1SShri Abhyankar #undef __FUNCT__
12975f2c45f1SShri Abhyankar #define __FUNCT__ "DMLocalToGlobalBegin_Network"
12985f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
12995f2c45f1SShri Abhyankar {
13005f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13015f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
13025f2c45f1SShri Abhyankar 
13035f2c45f1SShri Abhyankar   PetscFunctionBegin;
13045f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
13055f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13065f2c45f1SShri Abhyankar }
13075f2c45f1SShri Abhyankar 
13085f2c45f1SShri Abhyankar #undef __FUNCT__
13095f2c45f1SShri Abhyankar #define __FUNCT__ "DMLocalToGlobalEnd_Network"
13105f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
13115f2c45f1SShri Abhyankar {
13125f2c45f1SShri Abhyankar   PetscErrorCode ierr;
13135f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
13145f2c45f1SShri Abhyankar 
13155f2c45f1SShri Abhyankar   PetscFunctionBegin;
13165f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
13175f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
13185f2c45f1SShri Abhyankar }
1319