xref: /petsc/src/dm/impls/network/network.c (revision 7b6afd5b8febdfe9eecd254e43b0d168d52e271f)
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 /*@
6556ed216SShri Abhyankar   DMNetworkGetPlex - Gets the Plex DM associated with this network DM
7556ed216SShri Abhyankar 
8556ed216SShri Abhyankar   Not collective
9556ed216SShri Abhyankar 
10556ed216SShri Abhyankar   Input Parameters:
11556ed216SShri Abhyankar + netdm - the dm object
12556ed216SShri Abhyankar - plexmdm - the plex dm object
13556ed216SShri Abhyankar 
14556ed216SShri Abhyankar   Level: Advanced
15556ed216SShri Abhyankar 
16556ed216SShri Abhyankar .seealso: DMNetworkCreate()
17556ed216SShri Abhyankar @*/
18556ed216SShri Abhyankar PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm)
19556ed216SShri Abhyankar {
20556ed216SShri Abhyankar   DM_Network     *network = (DM_Network*) netdm->data;
21556ed216SShri Abhyankar 
22556ed216SShri Abhyankar   PetscFunctionBegin;
23556ed216SShri Abhyankar   *plexdm = network->plex;
24556ed216SShri Abhyankar   PetscFunctionReturn(0);
25556ed216SShri Abhyankar }
26556ed216SShri Abhyankar 
27556ed216SShri Abhyankar /*@
285f2c45f1SShri Abhyankar   DMNetworkSetSizes - Sets the local and global vertices and edges.
295f2c45f1SShri Abhyankar 
305f2c45f1SShri Abhyankar   Collective on DM
315f2c45f1SShri Abhyankar 
325f2c45f1SShri Abhyankar   Input Parameters:
335f2c45f1SShri Abhyankar + dm - the dm object
345f2c45f1SShri Abhyankar . nV - number of local vertices
355f2c45f1SShri Abhyankar . nE - number of local edges
365f2c45f1SShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE)
375f2c45f1SShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE)
385f2c45f1SShri Abhyankar 
395f2c45f1SShri Abhyankar    Notes
405f2c45f1SShri Abhyankar    If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang.
415f2c45f1SShri Abhyankar 
425f2c45f1SShri Abhyankar    You cannot change the sizes once they have been set
435f2c45f1SShri Abhyankar 
441b266c99SBarry Smith    Level: intermediate
451b266c99SBarry Smith 
461b266c99SBarry Smith .seealso: DMNetworkCreate()
475f2c45f1SShri Abhyankar @*/
485f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE)
495f2c45f1SShri Abhyankar {
505f2c45f1SShri Abhyankar   PetscErrorCode ierr;
515f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
525f2c45f1SShri Abhyankar   PetscInt       a[2],b[2];
535f2c45f1SShri Abhyankar 
545f2c45f1SShri Abhyankar   PetscFunctionBegin;
555f2c45f1SShri Abhyankar   PetscValidHeaderSpecific(dm,DM_CLASSID,1);
565f2c45f1SShri Abhyankar   if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4);
575f2c45f1SShri Abhyankar   if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5);
585f2c45f1SShri 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);
595f2c45f1SShri 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);
605f2c45f1SShri 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);
615f2c45f1SShri 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);
625f2c45f1SShri Abhyankar   if (NE < 0 || NV < 0) {
635f2c45f1SShri Abhyankar     a[0] = nV; a[1] = nE;
64b2566f29SBarry Smith     ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr);
655f2c45f1SShri Abhyankar     NV = b[0]; NE = b[1];
665f2c45f1SShri Abhyankar   }
675f2c45f1SShri Abhyankar   network->nNodes = nV;
685f2c45f1SShri Abhyankar   network->NNodes = NV;
695f2c45f1SShri Abhyankar   network->nEdges = nE;
705f2c45f1SShri Abhyankar   network->NEdges = NE;
715f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
725f2c45f1SShri Abhyankar }
735f2c45f1SShri Abhyankar 
745f2c45f1SShri Abhyankar /*@
755f2c45f1SShri Abhyankar   DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network
765f2c45f1SShri Abhyankar 
775f2c45f1SShri Abhyankar   Logically collective on DM
785f2c45f1SShri Abhyankar 
795f2c45f1SShri Abhyankar   Input Parameters:
805f2c45f1SShri Abhyankar . edges - list of edges
815f2c45f1SShri Abhyankar 
825f2c45f1SShri Abhyankar   Notes:
835f2c45f1SShri Abhyankar   There is no copy involved in this operation, only the pointer is referenced. The edgelist should
845f2c45f1SShri Abhyankar   not be destroyed before the call to DMNetworkLayoutSetUp
855f2c45f1SShri Abhyankar 
865f2c45f1SShri Abhyankar   Level: intermediate
875f2c45f1SShri Abhyankar 
885f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes
895f2c45f1SShri Abhyankar @*/
905f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[])
915f2c45f1SShri Abhyankar {
925f2c45f1SShri Abhyankar   DM_Network *network = (DM_Network*) dm->data;
935f2c45f1SShri Abhyankar 
945f2c45f1SShri Abhyankar   PetscFunctionBegin;
955f2c45f1SShri Abhyankar   network->edges = edgelist;
965f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
975f2c45f1SShri Abhyankar }
985f2c45f1SShri Abhyankar 
995f2c45f1SShri Abhyankar /*@
1005f2c45f1SShri Abhyankar   DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network
1015f2c45f1SShri Abhyankar 
1025f2c45f1SShri Abhyankar   Collective on DM
1035f2c45f1SShri Abhyankar 
1045f2c45f1SShri Abhyankar   Input Parameters
1055f2c45f1SShri Abhyankar . DM - the dmnetwork object
1065f2c45f1SShri Abhyankar 
1075f2c45f1SShri Abhyankar   Notes:
1085f2c45f1SShri Abhyankar   This routine should be called after the network sizes and edgelists have been provided. It creates
1095f2c45f1SShri Abhyankar   the bare layout of the network and sets up the network to begin insertion of components.
1105f2c45f1SShri Abhyankar 
1115f2c45f1SShri Abhyankar   All the components should be registered before calling this routine.
1125f2c45f1SShri Abhyankar 
1135f2c45f1SShri Abhyankar   Level: intermediate
1145f2c45f1SShri Abhyankar 
1155f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList
1165f2c45f1SShri Abhyankar @*/
1175f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm)
1185f2c45f1SShri Abhyankar {
1195f2c45f1SShri Abhyankar   PetscErrorCode ierr;
1205f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
1215f2c45f1SShri Abhyankar   PetscInt       dim = 1; /* One dimensional network */
1225f2c45f1SShri Abhyankar   PetscInt       numCorners=2;
1235f2c45f1SShri Abhyankar   PetscInt       spacedim=2;
1245f2c45f1SShri Abhyankar   double         *vertexcoords=NULL;
1255f2c45f1SShri Abhyankar   PetscInt       i;
1265f2c45f1SShri Abhyankar   PetscInt       ndata;
1275f2c45f1SShri Abhyankar 
1285f2c45f1SShri Abhyankar   PetscFunctionBegin;
1295f2c45f1SShri Abhyankar   if (network->nNodes) {
1309045477aSSatish Balay     ierr = PetscCalloc1(numCorners*network->nNodes,&vertexcoords);CHKERRQ(ierr);
1315f2c45f1SShri Abhyankar   }
1325f2c45f1SShri Abhyankar   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
1335f2c45f1SShri Abhyankar   if (network->nNodes) {
1345f2c45f1SShri Abhyankar     ierr = PetscFree(vertexcoords);CHKERRQ(ierr);
1355f2c45f1SShri Abhyankar   }
1365f2c45f1SShri Abhyankar   ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr);
1375f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr);
1385f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr);
1395f2c45f1SShri Abhyankar 
1405f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr);
1415f2c45f1SShri Abhyankar   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr);
1425f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr);
1435f2c45f1SShri Abhyankar   ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr);
1445f2c45f1SShri Abhyankar 
14524121865SAdrian Maldonado 
14624121865SAdrian Maldonado 
1475f2c45f1SShri Abhyankar   network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
1486caa05f4SBarry Smith   ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr);
1495f2c45f1SShri Abhyankar   for (i = network->pStart; i < network->pEnd; i++) {
150*7b6afd5bSHong Zhang     if (i < network->vStart) {
15162aa33baSHong Zhang       network->header[i].id = i;
152*7b6afd5bSHong Zhang     } else {
153*7b6afd5bSHong Zhang       network->header[i].id = i - network->vStart;
154*7b6afd5bSHong Zhang     }
1555f2c45f1SShri Abhyankar     network->header[i].ndata = 0;
1565f2c45f1SShri Abhyankar     ndata = network->header[i].ndata;
1575f2c45f1SShri Abhyankar     ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr);
1585f2c45f1SShri Abhyankar     network->header[i].offset[ndata] = 0;
1595f2c45f1SShri Abhyankar   }
160854ce69bSBarry Smith   ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr);
1615f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
1625f2c45f1SShri Abhyankar }
1635f2c45f1SShri Abhyankar 
16494ef8ddeSSatish Balay /*@C
1655f2c45f1SShri Abhyankar   DMNetworkRegisterComponent - Registers the network component
1665f2c45f1SShri Abhyankar 
1675f2c45f1SShri Abhyankar   Logically collective on DM
1685f2c45f1SShri Abhyankar 
1695f2c45f1SShri Abhyankar   Input Parameters
1705f2c45f1SShri Abhyankar + dm   - the network object
1715f2c45f1SShri Abhyankar . name - the component name
1725f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data
1735f2c45f1SShri Abhyankar 
1745f2c45f1SShri Abhyankar    Output Parameters
1755f2c45f1SShri Abhyankar .   key - an integer key that defines the component
1765f2c45f1SShri Abhyankar 
1775f2c45f1SShri Abhyankar    Notes
1785f2c45f1SShri Abhyankar    This routine should be called by all processors before calling DMNetworkLayoutSetup().
1795f2c45f1SShri Abhyankar 
1805f2c45f1SShri Abhyankar    Level: intermediate
1815f2c45f1SShri Abhyankar 
1825f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate
1835f2c45f1SShri Abhyankar @*/
1845f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key)
1855f2c45f1SShri Abhyankar {
1865f2c45f1SShri Abhyankar   PetscErrorCode        ierr;
1875f2c45f1SShri Abhyankar   DM_Network            *network = (DM_Network*) dm->data;
1885f2c45f1SShri Abhyankar   DMNetworkComponent    *component=&network->component[network->ncomponent];
1895f2c45f1SShri Abhyankar   PetscBool             flg=PETSC_FALSE;
1905f2c45f1SShri Abhyankar   PetscInt              i;
1915f2c45f1SShri Abhyankar 
1925f2c45f1SShri Abhyankar   PetscFunctionBegin;
1935f2c45f1SShri Abhyankar 
1945f2c45f1SShri Abhyankar   for (i=0; i < network->ncomponent; i++) {
1955f2c45f1SShri Abhyankar     ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr);
1965f2c45f1SShri Abhyankar     if (flg) {
1975f2c45f1SShri Abhyankar       *key = i;
1985f2c45f1SShri Abhyankar       PetscFunctionReturn(0);
1995f2c45f1SShri Abhyankar     }
2005f2c45f1SShri Abhyankar   }
2015f2c45f1SShri Abhyankar 
2025f2c45f1SShri Abhyankar   ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr);
2035f2c45f1SShri Abhyankar   component->size = size/sizeof(DMNetworkComponentGenericDataType);
2045f2c45f1SShri Abhyankar   *key = network->ncomponent;
2055f2c45f1SShri Abhyankar   network->ncomponent++;
2065f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2075f2c45f1SShri Abhyankar }
2085f2c45f1SShri Abhyankar 
2095f2c45f1SShri Abhyankar /*@
2105f2c45f1SShri Abhyankar   DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices.
2115f2c45f1SShri Abhyankar 
2125f2c45f1SShri Abhyankar   Not Collective
2135f2c45f1SShri Abhyankar 
2145f2c45f1SShri Abhyankar   Input Parameters:
2155f2c45f1SShri Abhyankar + dm - The DMNetwork object
2165f2c45f1SShri Abhyankar 
2175f2c45f1SShri Abhyankar   Output Paramters:
2185f2c45f1SShri Abhyankar + vStart - The first vertex point
2195f2c45f1SShri Abhyankar - vEnd   - One beyond the last vertex point
2205f2c45f1SShri Abhyankar 
2215f2c45f1SShri Abhyankar   Level: intermediate
2225f2c45f1SShri Abhyankar 
2235f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange
2245f2c45f1SShri Abhyankar @*/
2255f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd)
2265f2c45f1SShri Abhyankar {
2275f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
2285f2c45f1SShri Abhyankar 
2295f2c45f1SShri Abhyankar   PetscFunctionBegin;
2305f2c45f1SShri Abhyankar   if (vStart) *vStart = network->vStart;
2315f2c45f1SShri Abhyankar   if (vEnd) *vEnd = network->vEnd;
2325f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2335f2c45f1SShri Abhyankar }
2345f2c45f1SShri Abhyankar 
2355f2c45f1SShri Abhyankar /*@
2365f2c45f1SShri Abhyankar   DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges.
2375f2c45f1SShri Abhyankar 
2385f2c45f1SShri Abhyankar   Not Collective
2395f2c45f1SShri Abhyankar 
2405f2c45f1SShri Abhyankar   Input Parameters:
2415f2c45f1SShri Abhyankar + dm - The DMNetwork object
2425f2c45f1SShri Abhyankar 
2435f2c45f1SShri Abhyankar   Output Paramters:
2445f2c45f1SShri Abhyankar + eStart - The first edge point
2455f2c45f1SShri Abhyankar - eEnd   - One beyond the last edge point
2465f2c45f1SShri Abhyankar 
2475f2c45f1SShri Abhyankar   Level: intermediate
2485f2c45f1SShri Abhyankar 
2495f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange
2505f2c45f1SShri Abhyankar @*/
2515f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd)
2525f2c45f1SShri Abhyankar {
2535f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
2545f2c45f1SShri Abhyankar 
2555f2c45f1SShri Abhyankar   PetscFunctionBegin;
2565f2c45f1SShri Abhyankar   if (eStart) *eStart = network->eStart;
2575f2c45f1SShri Abhyankar   if (eEnd) *eEnd = network->eEnd;
2585f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
2595f2c45f1SShri Abhyankar }
2605f2c45f1SShri Abhyankar 
261*7b6afd5bSHong Zhang 
262*7b6afd5bSHong Zhang /*@
263*7b6afd5bSHong Zhang   DMNetworkGetEdgeVertexID - Get the user numbering for the edge or vertex.
264*7b6afd5bSHong Zhang 
265*7b6afd5bSHong Zhang   Not Collective
266*7b6afd5bSHong Zhang 
267*7b6afd5bSHong Zhang   Input Parameters:
268*7b6afd5bSHong Zhang + dm - DMNetwork object
269*7b6afd5bSHong Zhang - p  - vertex/edge point
270*7b6afd5bSHong Zhang 
271*7b6afd5bSHong Zhang   Output Paramters:
272*7b6afd5bSHong Zhang . id - user numbering
273*7b6afd5bSHong Zhang 
274*7b6afd5bSHong Zhang   Level: intermediate
275*7b6afd5bSHong Zhang 
276*7b6afd5bSHong Zhang .seealso: DMNetworkGetEdgeRange, DMNetworkGetVertexRange
277*7b6afd5bSHong Zhang @*/
278*7b6afd5bSHong Zhang PetscErrorCode DMNetworkGetEdgeVertexID(DM dm,PetscInt p,PetscInt *id)
279*7b6afd5bSHong Zhang {
280*7b6afd5bSHong Zhang   PetscErrorCode    ierr;
281*7b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
282*7b6afd5bSHong Zhang   PetscInt          offsetp;
283*7b6afd5bSHong Zhang   DMNetworkComponentHeader header;
284*7b6afd5bSHong Zhang 
285*7b6afd5bSHong Zhang   PetscFunctionBegin;
286*7b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
287*7b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
288*7b6afd5bSHong Zhang   *id    = header->id;
289*7b6afd5bSHong Zhang   PetscFunctionReturn(0);
290*7b6afd5bSHong Zhang }
291*7b6afd5bSHong Zhang 
2925f2c45f1SShri Abhyankar /*@
2935f2c45f1SShri Abhyankar   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
2945f2c45f1SShri Abhyankar 
2955f2c45f1SShri Abhyankar   Not Collective
2965f2c45f1SShri Abhyankar 
2975f2c45f1SShri Abhyankar   Input Parameters:
2985f2c45f1SShri Abhyankar + dm           - The DMNetwork object
2995f2c45f1SShri Abhyankar . p            - vertex/edge point
3005f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
3015f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
3025f2c45f1SShri Abhyankar 
3035f2c45f1SShri Abhyankar   Level: intermediate
3045f2c45f1SShri Abhyankar 
3055f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
3065f2c45f1SShri Abhyankar @*/
3075f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
3085f2c45f1SShri Abhyankar {
3095f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
31043a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
3115f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
3125f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
3135f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
3145f2c45f1SShri Abhyankar 
3155f2c45f1SShri Abhyankar   PetscFunctionBegin;
316fa58f0a9SHong 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);
317fa58f0a9SHong Zhang 
31843a39a44SBarry Smith   header->size[header->ndata] = component->size;
31943a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
3205f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
3215f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
3225f2c45f1SShri Abhyankar 
3235f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
3245f2c45f1SShri Abhyankar   header->ndata++;
3255f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3265f2c45f1SShri Abhyankar }
3275f2c45f1SShri Abhyankar 
3285f2c45f1SShri Abhyankar /*@
3295f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
3305f2c45f1SShri Abhyankar 
3315f2c45f1SShri Abhyankar   Not Collective
3325f2c45f1SShri Abhyankar 
3335f2c45f1SShri Abhyankar   Input Parameters:
3345f2c45f1SShri Abhyankar + dm - The DMNetwork object
3355f2c45f1SShri Abhyankar . p  - vertex/edge point
3365f2c45f1SShri Abhyankar 
3375f2c45f1SShri Abhyankar   Output Parameters:
3385f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
3395f2c45f1SShri Abhyankar 
3405f2c45f1SShri Abhyankar   Level: intermediate
3415f2c45f1SShri Abhyankar 
3425f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
3435f2c45f1SShri Abhyankar @*/
3445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
3455f2c45f1SShri Abhyankar {
3465f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3475f2c45f1SShri Abhyankar   PetscInt       offset;
3485f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3495f2c45f1SShri Abhyankar 
3505f2c45f1SShri Abhyankar   PetscFunctionBegin;
3515f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
3525f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
3535f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3545f2c45f1SShri Abhyankar }
3555f2c45f1SShri Abhyankar 
3565f2c45f1SShri Abhyankar /*@
3575f2c45f1SShri Abhyankar   DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the
3585f2c45f1SShri Abhyankar                                     component value from the component data array
3595f2c45f1SShri Abhyankar 
3605f2c45f1SShri Abhyankar   Not Collective
3615f2c45f1SShri Abhyankar 
3625f2c45f1SShri Abhyankar   Input Parameters:
3635f2c45f1SShri Abhyankar + dm      - The DMNetwork object
3645f2c45f1SShri Abhyankar . p       - vertex/edge point
3655f2c45f1SShri Abhyankar - compnum - component number
3665f2c45f1SShri Abhyankar 
3675f2c45f1SShri Abhyankar   Output Parameters:
3685f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component
3695f2c45f1SShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
3705f2c45f1SShri Abhyankar 
3715f2c45f1SShri Abhyankar   Notes:
3725f2c45f1SShri Abhyankar   Typical usage:
3735f2c45f1SShri Abhyankar 
3745f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
3755f2c45f1SShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
3765f2c45f1SShri Abhyankar   Loop over vertices or edges
3775f2c45f1SShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
3785f2c45f1SShri Abhyankar     Loop over numcomps
3795f2c45f1SShri Abhyankar       DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset);
3805f2c45f1SShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
3815f2c45f1SShri Abhyankar 
3825f2c45f1SShri Abhyankar   Level: intermediate
3835f2c45f1SShri Abhyankar 
3845f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
3855f2c45f1SShri Abhyankar @*/
3865f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
3875f2c45f1SShri Abhyankar {
3885f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
3895f2c45f1SShri Abhyankar   PetscInt                 offsetp;
3905f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
3915f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
3925f2c45f1SShri Abhyankar 
3935f2c45f1SShri Abhyankar   PetscFunctionBegin;
3945f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
3955f2c45f1SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
396d36f4e81SHong Zhang   if (compkey) *compkey = header->key[compnum];
397d36f4e81SHong Zhang   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
3985f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3995f2c45f1SShri Abhyankar }
4005f2c45f1SShri Abhyankar 
4015f2c45f1SShri Abhyankar /*@
4025f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
4035f2c45f1SShri Abhyankar 
4045f2c45f1SShri Abhyankar   Not Collective
4055f2c45f1SShri Abhyankar 
4065f2c45f1SShri Abhyankar   Input Parameters:
4075f2c45f1SShri Abhyankar + dm     - The DMNetwork object
4085f2c45f1SShri Abhyankar - p      - the edge/vertex point
4095f2c45f1SShri Abhyankar 
4105f2c45f1SShri Abhyankar   Output Parameters:
4115f2c45f1SShri Abhyankar . offset - the offset
4125f2c45f1SShri Abhyankar 
4135f2c45f1SShri Abhyankar   Level: intermediate
4145f2c45f1SShri Abhyankar 
4155f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
4165f2c45f1SShri Abhyankar @*/
4175f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
4185f2c45f1SShri Abhyankar {
4195f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4205f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4215f2c45f1SShri Abhyankar 
4225f2c45f1SShri Abhyankar   PetscFunctionBegin;
4235f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
4245f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4255f2c45f1SShri Abhyankar }
4265f2c45f1SShri Abhyankar 
4275f2c45f1SShri Abhyankar /*@
4285f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
4295f2c45f1SShri Abhyankar 
4305f2c45f1SShri Abhyankar   Not Collective
4315f2c45f1SShri Abhyankar 
4325f2c45f1SShri Abhyankar   Input Parameters:
4335f2c45f1SShri Abhyankar + dm      - The DMNetwork object
4345f2c45f1SShri Abhyankar - p       - the edge/vertex point
4355f2c45f1SShri Abhyankar 
4365f2c45f1SShri Abhyankar   Output Parameters:
4375f2c45f1SShri Abhyankar . offsetg - the offset
4385f2c45f1SShri Abhyankar 
4395f2c45f1SShri Abhyankar   Level: intermediate
4405f2c45f1SShri Abhyankar 
4415f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
4425f2c45f1SShri Abhyankar @*/
4435f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
4445f2c45f1SShri Abhyankar {
4455f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4465f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4475f2c45f1SShri Abhyankar 
4485f2c45f1SShri Abhyankar   PetscFunctionBegin;
4495f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
450dd6f46cdSHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost node */
4515f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4525f2c45f1SShri Abhyankar }
4535f2c45f1SShri Abhyankar 
45424121865SAdrian Maldonado /*@
45524121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
45624121865SAdrian Maldonado 
45724121865SAdrian Maldonado   Not Collective
45824121865SAdrian Maldonado 
45924121865SAdrian Maldonado   Input Parameters:
46024121865SAdrian Maldonado + dm     - The DMNetwork object
46124121865SAdrian Maldonado - p      - the edge point
46224121865SAdrian Maldonado 
46324121865SAdrian Maldonado   Output Parameters:
46424121865SAdrian Maldonado . offset - the offset
46524121865SAdrian Maldonado 
46624121865SAdrian Maldonado   Level: intermediate
46724121865SAdrian Maldonado 
46824121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
46924121865SAdrian Maldonado @*/
47024121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
47124121865SAdrian Maldonado {
47224121865SAdrian Maldonado   PetscErrorCode ierr;
47324121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
47424121865SAdrian Maldonado 
47524121865SAdrian Maldonado   PetscFunctionBegin;
47624121865SAdrian Maldonado 
47724121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
47824121865SAdrian Maldonado   PetscFunctionReturn(0);
47924121865SAdrian Maldonado }
48024121865SAdrian Maldonado 
48124121865SAdrian Maldonado /*@
48224121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
48324121865SAdrian Maldonado 
48424121865SAdrian Maldonado   Not Collective
48524121865SAdrian Maldonado 
48624121865SAdrian Maldonado   Input Parameters:
48724121865SAdrian Maldonado + dm     - The DMNetwork object
48824121865SAdrian Maldonado - p      - the vertex point
48924121865SAdrian Maldonado 
49024121865SAdrian Maldonado   Output Parameters:
49124121865SAdrian Maldonado . offset - the offset
49224121865SAdrian Maldonado 
49324121865SAdrian Maldonado   Level: intermediate
49424121865SAdrian Maldonado 
49524121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
49624121865SAdrian Maldonado @*/
49724121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
49824121865SAdrian Maldonado {
49924121865SAdrian Maldonado   PetscErrorCode ierr;
50024121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
50124121865SAdrian Maldonado 
50224121865SAdrian Maldonado   PetscFunctionBegin;
50324121865SAdrian Maldonado 
50424121865SAdrian Maldonado   p -= network->vStart;
50524121865SAdrian Maldonado 
50624121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
50724121865SAdrian Maldonado   PetscFunctionReturn(0);
50824121865SAdrian Maldonado }
5095f2c45f1SShri Abhyankar /*@
5105f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
5115f2c45f1SShri Abhyankar 
5125f2c45f1SShri Abhyankar   Not Collective
5135f2c45f1SShri Abhyankar 
5145f2c45f1SShri Abhyankar   Input Parameters:
5155f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
5165f2c45f1SShri Abhyankar . p    - the vertex/edge point
5175f2c45f1SShri Abhyankar - nvar - number of additional variables
5185f2c45f1SShri Abhyankar 
5195f2c45f1SShri Abhyankar   Level: intermediate
5205f2c45f1SShri Abhyankar 
5215f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
5225f2c45f1SShri Abhyankar @*/
5235f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
5245f2c45f1SShri Abhyankar {
5255f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5265f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5275f2c45f1SShri Abhyankar 
5285f2c45f1SShri Abhyankar   PetscFunctionBegin;
5295f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
5305f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5315f2c45f1SShri Abhyankar }
5325f2c45f1SShri Abhyankar 
53327f51fceSHong Zhang /*@
53427f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
53527f51fceSHong Zhang 
53627f51fceSHong Zhang   Not Collective
53727f51fceSHong Zhang 
53827f51fceSHong Zhang   Input Parameters:
53927f51fceSHong Zhang + dm   - The DMNetworkObject
54027f51fceSHong Zhang - p    - the vertex/edge point
54127f51fceSHong Zhang 
54227f51fceSHong Zhang   Output Parameters:
54327f51fceSHong Zhang . nvar - number of variables
54427f51fceSHong Zhang 
54527f51fceSHong Zhang   Level: intermediate
54627f51fceSHong Zhang 
54727f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
54827f51fceSHong Zhang @*/
54927f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
55027f51fceSHong Zhang {
55127f51fceSHong Zhang   PetscErrorCode ierr;
55227f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
55327f51fceSHong Zhang 
55427f51fceSHong Zhang   PetscFunctionBegin;
55527f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
55627f51fceSHong Zhang   PetscFunctionReturn(0);
55727f51fceSHong Zhang }
55827f51fceSHong Zhang 
5595f2c45f1SShri Abhyankar /*@
5605f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
5615f2c45f1SShri Abhyankar 
5625f2c45f1SShri Abhyankar   Not Collective
5635f2c45f1SShri Abhyankar 
5645f2c45f1SShri Abhyankar   Input Parameters:
5655f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
5665f2c45f1SShri Abhyankar . p    - the vertex/edge point
5675f2c45f1SShri Abhyankar - nvar - number of variables
5685f2c45f1SShri Abhyankar 
5695f2c45f1SShri Abhyankar   Level: intermediate
5705f2c45f1SShri Abhyankar 
5715f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
5725f2c45f1SShri Abhyankar @*/
5735f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
5745f2c45f1SShri Abhyankar {
5755f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5765f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5775f2c45f1SShri Abhyankar 
5785f2c45f1SShri Abhyankar   PetscFunctionBegin;
5795f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
5805f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5815f2c45f1SShri Abhyankar }
5825f2c45f1SShri Abhyankar 
5835f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
5845f2c45f1SShri Abhyankar    function is called during DMSetUp() */
5855f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
5865f2c45f1SShri Abhyankar {
5875f2c45f1SShri Abhyankar   PetscErrorCode              ierr;
5885f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5895f2c45f1SShri Abhyankar   PetscInt                    arr_size;
5905f2c45f1SShri Abhyankar   PetscInt                    p,offset,offsetp;
5915f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
5925f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
5935f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType      *componentdataarray;
5945f2c45f1SShri Abhyankar   PetscInt ncomp, i;
5955f2c45f1SShri Abhyankar 
5965f2c45f1SShri Abhyankar   PetscFunctionBegin;
5975f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
5985f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
59975b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
6005f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
6015f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
6025f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
6035f2c45f1SShri Abhyankar     /* Copy header */
6045f2c45f1SShri Abhyankar     header = &network->header[p];
605302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
6065f2c45f1SShri Abhyankar     /* Copy data */
6075f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
6085f2c45f1SShri Abhyankar     ncomp = header->ndata;
6095f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
6105f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
611302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
6125f2c45f1SShri Abhyankar     }
6135f2c45f1SShri Abhyankar   }
6145f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6155f2c45f1SShri Abhyankar }
6165f2c45f1SShri Abhyankar 
6175f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
6185f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
6195f2c45f1SShri Abhyankar {
6205f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6215f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6225f2c45f1SShri Abhyankar 
6235f2c45f1SShri Abhyankar   PetscFunctionBegin;
6245f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
6255f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6265f2c45f1SShri Abhyankar }
6275f2c45f1SShri Abhyankar 
6285f2c45f1SShri Abhyankar /*@C
6295f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
6305f2c45f1SShri Abhyankar 
6315f2c45f1SShri Abhyankar   Not Collective
6325f2c45f1SShri Abhyankar 
6335f2c45f1SShri Abhyankar   Input Parameters:
6345f2c45f1SShri Abhyankar . dm - The DMNetwork Object
6355f2c45f1SShri Abhyankar 
6365f2c45f1SShri Abhyankar   Output Parameters:
6375f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
6385f2c45f1SShri Abhyankar 
6395f2c45f1SShri Abhyankar   Level: intermediate
6405f2c45f1SShri Abhyankar 
6415f2c45f1SShri Abhyankar .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents
6425f2c45f1SShri Abhyankar @*/
6435f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
6445f2c45f1SShri Abhyankar {
6455f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6465f2c45f1SShri Abhyankar 
6475f2c45f1SShri Abhyankar   PetscFunctionBegin;
6485f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
6495f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6505f2c45f1SShri Abhyankar }
6515f2c45f1SShri Abhyankar 
65224121865SAdrian Maldonado /* Get a subsection from a range of points */
65324121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
65424121865SAdrian Maldonado {
65524121865SAdrian Maldonado   PetscErrorCode ierr;
65624121865SAdrian Maldonado   PetscInt       i, nvar;
65724121865SAdrian Maldonado 
65824121865SAdrian Maldonado   PetscFunctionBegin;
65924121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
66024121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
66124121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
66224121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
66324121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
66424121865SAdrian Maldonado   }
66524121865SAdrian Maldonado 
66624121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
66724121865SAdrian Maldonado   PetscFunctionReturn(0);
66824121865SAdrian Maldonado }
66924121865SAdrian Maldonado 
67024121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
67124121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
67224121865SAdrian Maldonado {
67324121865SAdrian Maldonado   PetscErrorCode ierr;
67424121865SAdrian Maldonado   PetscInt       i, *subpoints;
67524121865SAdrian Maldonado 
67624121865SAdrian Maldonado   PetscFunctionBegin;
67724121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
67824121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
67924121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
68024121865SAdrian Maldonado     subpoints[i - pstart] = i;
68124121865SAdrian Maldonado   }
682459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
68324121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
68424121865SAdrian Maldonado   PetscFunctionReturn(0);
68524121865SAdrian Maldonado }
68624121865SAdrian Maldonado 
68724121865SAdrian Maldonado /*@
68824121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
68924121865SAdrian Maldonado 
69024121865SAdrian Maldonado   Collective
69124121865SAdrian Maldonado 
69224121865SAdrian Maldonado   Input Parameters:
69324121865SAdrian Maldonado . dm   - The DMNetworkObject
69424121865SAdrian Maldonado 
69524121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
69624121865SAdrian Maldonado 
69724121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
69824121865SAdrian Maldonado 
69924121865SAdrian Maldonado   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
70024121865SAdrian Maldonado 
70124121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
70224121865SAdrian Maldonado 
70324121865SAdrian Maldonado   Level: intermediate
70424121865SAdrian Maldonado 
70524121865SAdrian Maldonado @*/
70624121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
70724121865SAdrian Maldonado {
70824121865SAdrian Maldonado   PetscErrorCode ierr;
70924121865SAdrian Maldonado   MPI_Comm       comm;
7109852e123SBarry Smith   PetscMPIInt    rank, size;
71124121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
71224121865SAdrian Maldonado 
713eab1376dSHong Zhang   PetscFunctionBegin;
71424121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
71524121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
7169852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
71724121865SAdrian Maldonado 
71824121865SAdrian Maldonado   /* Create maps for vertices and edges */
71924121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
72024121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
72124121865SAdrian Maldonado 
72224121865SAdrian Maldonado   /* Create local sub-sections */
72324121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
72424121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
72524121865SAdrian Maldonado 
7269852e123SBarry Smith   if (size > 1) {
72724121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
72824121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
72924121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
73024121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
73124121865SAdrian Maldonado   } else {
73224121865SAdrian Maldonado   /* create structures for vertex */
73324121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
73424121865SAdrian Maldonado   /* create structures for edge */
73524121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
73624121865SAdrian Maldonado   }
73724121865SAdrian Maldonado 
73824121865SAdrian Maldonado 
73924121865SAdrian Maldonado   /* Add viewers */
74024121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
74124121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
74224121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
74324121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
74424121865SAdrian Maldonado 
74524121865SAdrian Maldonado   PetscFunctionReturn(0);
74624121865SAdrian Maldonado }
747*7b6afd5bSHong Zhang 
7485f2c45f1SShri Abhyankar /*@
7495f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
7505f2c45f1SShri Abhyankar 
7515f2c45f1SShri Abhyankar   Collective
7525f2c45f1SShri Abhyankar 
7535f2c45f1SShri Abhyankar   Input Parameter:
754d3464fd4SAdrian Maldonado + DM - the DMNetwork object
7555f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
7565f2c45f1SShri Abhyankar 
7575f2c45f1SShri Abhyankar   Notes:
7588b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
7595f2c45f1SShri Abhyankar 
7605f2c45f1SShri Abhyankar   Level: intermediate
7615f2c45f1SShri Abhyankar 
7625f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
7635f2c45f1SShri Abhyankar @*/
764d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
7655f2c45f1SShri Abhyankar {
766d3464fd4SAdrian Maldonado   MPI_Comm       comm;
7675f2c45f1SShri Abhyankar   PetscErrorCode ierr;
768d3464fd4SAdrian Maldonado   PetscMPIInt    size;
769d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
770d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
7715f2c45f1SShri Abhyankar   PetscSF        pointsf;
7725f2c45f1SShri Abhyankar   DM             newDM;
77351ac5effSHong Zhang   PetscPartitioner part;
7745f2c45f1SShri Abhyankar 
7755f2c45f1SShri Abhyankar   PetscFunctionBegin;
776d3464fd4SAdrian Maldonado 
777d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
778d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
779d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
780d3464fd4SAdrian Maldonado 
781d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
7825f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
7835f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
78451ac5effSHong Zhang 
78551ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
78651ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
78751ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
78851ac5effSHong Zhang 
7895f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
79080cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
79151ac5effSHong Zhang 
7925f2c45f1SShri Abhyankar   /* Distribute dof section */
793d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
7945f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
795d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
79651ac5effSHong Zhang 
7975f2c45f1SShri Abhyankar   /* Distribute data and associated section */
79831da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
79924121865SAdrian Maldonado 
8005f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
8015f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
8025f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
8035f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
8045f2c45f1SShri Abhyankar   newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart;
8055f2c45f1SShri Abhyankar   newDMnetwork->NNodes = oldDMnetwork->NNodes;
8065f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
80724121865SAdrian Maldonado 
8085f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
8095f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
8105f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
8115f2c45f1SShri Abhyankar 
81224121865SAdrian Maldonado   /* Destroy point SF */
81324121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
81424121865SAdrian Maldonado 
815d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
816d3464fd4SAdrian Maldonado   *dm  = newDM;
8175f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8185f2c45f1SShri Abhyankar }
8195f2c45f1SShri Abhyankar 
82024121865SAdrian Maldonado /*@C
82124121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
82224121865SAdrian Maldonado 
82324121865SAdrian Maldonado   Input Parameters:
82424121865SAdrian Maldonado + masterSF - the original SF structure
82524121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
82624121865SAdrian Maldonado 
82724121865SAdrian Maldonado   Output Parameters:
82824121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
82924121865SAdrian Maldonado */
83024121865SAdrian Maldonado 
83124121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
83224121865SAdrian Maldonado 
83324121865SAdrian Maldonado   PetscErrorCode        ierr;
83424121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
83524121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
83624121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
83724121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
83824121865SAdrian Maldonado   const PetscInt        *ilocal;
83924121865SAdrian Maldonado   const PetscSFNode     *iremote;
84024121865SAdrian Maldonado 
84124121865SAdrian Maldonado   PetscFunctionBegin;
84224121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
84324121865SAdrian Maldonado 
84424121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
84524121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
84624121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
84724121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
84824121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
84924121865SAdrian Maldonado   }
85024121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
85124121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
85224121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
85324121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
85424121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
85524121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
85624121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
8574b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
8584b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
85924121865SAdrian Maldonado   nleaves_sub = 0;
86024121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
86124121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
86224121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
8634b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
86424121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
86524121865SAdrian Maldonado       nleaves_sub += 1;
86624121865SAdrian Maldonado     }
86724121865SAdrian Maldonado   }
86824121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
86924121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
87024121865SAdrian Maldonado 
87124121865SAdrian Maldonado   /* Create new subSF */
87224121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
87324121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
8744b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
87524121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
8764b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
87724121865SAdrian Maldonado 
87824121865SAdrian Maldonado   PetscFunctionReturn(0);
87924121865SAdrian Maldonado }
88024121865SAdrian Maldonado 
8815f2c45f1SShri Abhyankar /*@C
8825f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
8835f2c45f1SShri Abhyankar 
8845f2c45f1SShri Abhyankar   Not Collective
8855f2c45f1SShri Abhyankar 
8865f2c45f1SShri Abhyankar   Input Parameters:
8875f2c45f1SShri Abhyankar + dm - The DMNetwork object
8885f2c45f1SShri Abhyankar - p  - the vertex point
8895f2c45f1SShri Abhyankar 
8905f2c45f1SShri Abhyankar   Output Paramters:
8915f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
8925f2c45f1SShri Abhyankar - edges  - List of edge points
8935f2c45f1SShri Abhyankar 
8945f2c45f1SShri Abhyankar   Level: intermediate
8955f2c45f1SShri Abhyankar 
8965f2c45f1SShri Abhyankar   Fortran Notes:
8975f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
8985f2c45f1SShri Abhyankar   include petsc.h90 in your code.
8995f2c45f1SShri Abhyankar 
9005f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes
9015f2c45f1SShri Abhyankar @*/
9025f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
9035f2c45f1SShri Abhyankar {
9045f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9055f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9065f2c45f1SShri Abhyankar 
9075f2c45f1SShri Abhyankar   PetscFunctionBegin;
9085f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
9095f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
9105f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9115f2c45f1SShri Abhyankar }
9125f2c45f1SShri Abhyankar 
9135f2c45f1SShri Abhyankar /*@C
914596e729fSHong Zhang   DMNetworkGetConnectedNodes - Return the connected vertices for this edge point
9155f2c45f1SShri Abhyankar 
9165f2c45f1SShri Abhyankar   Not Collective
9175f2c45f1SShri Abhyankar 
9185f2c45f1SShri Abhyankar   Input Parameters:
9195f2c45f1SShri Abhyankar + dm - The DMNetwork object
9205f2c45f1SShri Abhyankar - p  - the edge point
9215f2c45f1SShri Abhyankar 
9225f2c45f1SShri Abhyankar   Output Paramters:
9235f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
9245f2c45f1SShri Abhyankar 
9255f2c45f1SShri Abhyankar   Level: intermediate
9265f2c45f1SShri Abhyankar 
9275f2c45f1SShri Abhyankar   Fortran Notes:
9285f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
9295f2c45f1SShri Abhyankar   include petsc.h90 in your code.
9305f2c45f1SShri Abhyankar 
9315f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
9325f2c45f1SShri Abhyankar @*/
9335f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[])
9345f2c45f1SShri Abhyankar {
9355f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9365f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9375f2c45f1SShri Abhyankar 
9385f2c45f1SShri Abhyankar   PetscFunctionBegin;
9395f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
9405f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9415f2c45f1SShri Abhyankar }
9425f2c45f1SShri Abhyankar 
9435f2c45f1SShri Abhyankar /*@
9445f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
9455f2c45f1SShri Abhyankar 
9465f2c45f1SShri Abhyankar   Not Collective
9475f2c45f1SShri Abhyankar 
9485f2c45f1SShri Abhyankar   Input Parameters:
9495f2c45f1SShri Abhyankar + dm - The DMNetwork object
9505f2c45f1SShri Abhyankar . p  - the vertex point
9515f2c45f1SShri Abhyankar 
9525f2c45f1SShri Abhyankar   Output Parameter:
9535f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
9545f2c45f1SShri Abhyankar 
9555f2c45f1SShri Abhyankar   Level: intermediate
9565f2c45f1SShri Abhyankar 
9575f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange
9585f2c45f1SShri Abhyankar @*/
9595f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
9605f2c45f1SShri Abhyankar {
9615f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9625f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9635f2c45f1SShri Abhyankar   PetscInt       offsetg;
9645f2c45f1SShri Abhyankar   PetscSection   sectiong;
9655f2c45f1SShri Abhyankar 
9665f2c45f1SShri Abhyankar   PetscFunctionBegin;
9675f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
9685f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
9695f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
9705f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
9715f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9725f2c45f1SShri Abhyankar }
9735f2c45f1SShri Abhyankar 
9745f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
9755f2c45f1SShri Abhyankar {
9765f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9775f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
9785f2c45f1SShri Abhyankar 
9795f2c45f1SShri Abhyankar   PetscFunctionBegin;
9805f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
9815f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
9825f2c45f1SShri Abhyankar 
9835f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
9845f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
9855f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9865f2c45f1SShri Abhyankar }
9875f2c45f1SShri Abhyankar 
9881ad426b7SHong Zhang /*@
98917df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
9901ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
9911ad426b7SHong Zhang 
9921ad426b7SHong Zhang     Collective
9931ad426b7SHong Zhang 
9941ad426b7SHong Zhang     Input Parameters:
99583b2e829SHong Zhang +   dm - The DMNetwork object
99683b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
99783b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
9981ad426b7SHong Zhang 
9991ad426b7SHong Zhang     Level: intermediate
10001ad426b7SHong Zhang 
10011ad426b7SHong Zhang @*/
100283b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
10031ad426b7SHong Zhang {
10041ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
10051ad426b7SHong Zhang 
10061ad426b7SHong Zhang   PetscFunctionBegin;
100783b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
100883b2e829SHong Zhang   network->userVertexJacobian = vflg;
10091ad426b7SHong Zhang   PetscFunctionReturn(0);
10101ad426b7SHong Zhang }
10111ad426b7SHong Zhang 
10121ad426b7SHong Zhang /*@
101383b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
101483b2e829SHong Zhang 
101583b2e829SHong Zhang     Not Collective
101683b2e829SHong Zhang 
101783b2e829SHong Zhang     Input Parameters:
101883b2e829SHong Zhang +   dm - The DMNetwork object
101983b2e829SHong Zhang .   p  - the edge point
10203e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
10213e97b6e8SHong Zhang         J[0]: this edge
10223e97b6e8SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes()
102383b2e829SHong Zhang 
102483b2e829SHong Zhang     Level: intermediate
102583b2e829SHong Zhang 
102683b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
102783b2e829SHong Zhang @*/
102883b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
102983b2e829SHong Zhang {
103083b2e829SHong Zhang   PetscErrorCode ierr;
103183b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
103283b2e829SHong Zhang 
103383b2e829SHong Zhang   PetscFunctionBegin;
1034883e35e8SHong Zhang   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
103583b2e829SHong Zhang   if (!network->Je) {
1036883e35e8SHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
103783b2e829SHong Zhang   }
1038883e35e8SHong Zhang   network->Je[3*p]   = J[0];
1039883e35e8SHong Zhang   network->Je[3*p+1] = J[1];
1040883e35e8SHong Zhang   network->Je[3*p+2] = J[2];
104183b2e829SHong Zhang   PetscFunctionReturn(0);
104283b2e829SHong Zhang }
104383b2e829SHong Zhang 
104483b2e829SHong Zhang /*@
104576ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
10461ad426b7SHong Zhang 
10471ad426b7SHong Zhang     Not Collective
10481ad426b7SHong Zhang 
10491ad426b7SHong Zhang     Input Parameters:
10501ad426b7SHong Zhang +   dm - The DMNetwork object
10511ad426b7SHong Zhang .   p  - the vertex point
10523e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
10533e97b6e8SHong Zhang         J[0]:       this vertex
10543e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
10553e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
10561ad426b7SHong Zhang 
10571ad426b7SHong Zhang     Level: intermediate
10581ad426b7SHong Zhang 
105983b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
10601ad426b7SHong Zhang @*/
1061883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
10625f2c45f1SShri Abhyankar {
10635f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10645f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
1065f4431b8cSHong Zhang   PetscInt       i,*vptr,nedges,vStart,vEnd;
1066883e35e8SHong Zhang   const PetscInt *edges;
10675f2c45f1SShri Abhyankar 
10685f2c45f1SShri Abhyankar   PetscFunctionBegin;
1069883e35e8SHong Zhang   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1070883e35e8SHong Zhang 
1071883e35e8SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1072f4431b8cSHong Zhang 
107383b2e829SHong Zhang   if (!network->Jv) {
1074f4431b8cSHong Zhang     PetscInt nNodes = network->nNodes,nedges_total;
1075883e35e8SHong Zhang     /* count nvertex_total */
10763e97b6e8SHong Zhang     nedges_total = 0;
1077883e35e8SHong Zhang     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1078f4431b8cSHong Zhang     ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr);
1079f4431b8cSHong Zhang 
1080883e35e8SHong Zhang     vptr[0] = 0;
1081f4431b8cSHong Zhang     for (i=0; i<nNodes; i++) {
1082f4431b8cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1083883e35e8SHong Zhang       nedges_total += nedges;
1084f4431b8cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
10851ad426b7SHong Zhang     }
10863e97b6e8SHong Zhang 
1087f4431b8cSHong Zhang     ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr);
1088883e35e8SHong Zhang     network->Jvptr = vptr;
1089883e35e8SHong Zhang   }
1090883e35e8SHong Zhang 
1091883e35e8SHong Zhang   vptr = network->Jvptr;
10923e97b6e8SHong Zhang   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
10933e97b6e8SHong Zhang 
10943e97b6e8SHong Zhang   /* Set Jacobian for each supporting edge and connected vertex */
1095883e35e8SHong Zhang   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1096883e35e8SHong Zhang   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1097883e35e8SHong Zhang   PetscFunctionReturn(0);
1098883e35e8SHong Zhang }
1099883e35e8SHong Zhang 
1100e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
11015cf7da58SHong Zhang {
11025cf7da58SHong Zhang   PetscErrorCode ierr;
11035cf7da58SHong Zhang   PetscInt       j;
11045cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
11055cf7da58SHong Zhang 
11065cf7da58SHong Zhang   PetscFunctionBegin;
11075cf7da58SHong Zhang   if (!ghost) {
11085cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11095cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11105cf7da58SHong Zhang     }
11115cf7da58SHong Zhang   } else {
11125cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11135cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11145cf7da58SHong Zhang     }
11155cf7da58SHong Zhang   }
11165cf7da58SHong Zhang   PetscFunctionReturn(0);
11175cf7da58SHong Zhang }
11185cf7da58SHong Zhang 
1119e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
11205cf7da58SHong Zhang {
11215cf7da58SHong Zhang   PetscErrorCode ierr;
11225cf7da58SHong Zhang   PetscInt       j,ncols_u;
11235cf7da58SHong Zhang   PetscScalar    val;
11245cf7da58SHong Zhang 
11255cf7da58SHong Zhang   PetscFunctionBegin;
11265cf7da58SHong Zhang   if (!ghost) {
11275cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11285cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
11295cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
11305cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11315cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
11325cf7da58SHong Zhang     }
11335cf7da58SHong Zhang   } else {
11345cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11355cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
11365cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
11375cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11385cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
11395cf7da58SHong Zhang     }
11405cf7da58SHong Zhang   }
11415cf7da58SHong Zhang   PetscFunctionReturn(0);
11425cf7da58SHong Zhang }
11435cf7da58SHong Zhang 
1144e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
11455cf7da58SHong Zhang {
11465cf7da58SHong Zhang   PetscErrorCode ierr;
11475cf7da58SHong Zhang 
11485cf7da58SHong Zhang   PetscFunctionBegin;
11495cf7da58SHong Zhang   if (Ju) {
11505cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
11515cf7da58SHong Zhang   } else {
11525cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
11535cf7da58SHong Zhang   }
11545cf7da58SHong Zhang   PetscFunctionReturn(0);
11555cf7da58SHong Zhang }
11565cf7da58SHong Zhang 
1157e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1158883e35e8SHong Zhang {
1159883e35e8SHong Zhang   PetscErrorCode ierr;
1160883e35e8SHong Zhang   PetscInt       j,*cols;
1161883e35e8SHong Zhang   PetscScalar    *zeros;
1162883e35e8SHong Zhang 
1163883e35e8SHong Zhang   PetscFunctionBegin;
1164883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1165883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1166883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1167883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
11681ad426b7SHong Zhang   PetscFunctionReturn(0);
11691ad426b7SHong Zhang }
1170a4e85ca8SHong Zhang 
1171e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
11723e97b6e8SHong Zhang {
11733e97b6e8SHong Zhang   PetscErrorCode ierr;
11743e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
11753e97b6e8SHong Zhang   const PetscInt *cols;
11763e97b6e8SHong Zhang   PetscScalar    zero=0.0;
11773e97b6e8SHong Zhang 
11783e97b6e8SHong Zhang   PetscFunctionBegin;
11793e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
11803e97b6e8SHong 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);
11813e97b6e8SHong Zhang 
11823e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
11833e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
11843e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
11853e97b6e8SHong Zhang       col = cols[j] + cstart;
11863e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
11873e97b6e8SHong Zhang     }
11883e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
11893e97b6e8SHong Zhang   }
11903e97b6e8SHong Zhang   PetscFunctionReturn(0);
11913e97b6e8SHong Zhang }
11921ad426b7SHong Zhang 
1193e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1194a4e85ca8SHong Zhang {
1195a4e85ca8SHong Zhang   PetscErrorCode ierr;
1196f4431b8cSHong Zhang 
1197a4e85ca8SHong Zhang   PetscFunctionBegin;
1198a4e85ca8SHong Zhang   if (Ju) {
1199a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1200a4e85ca8SHong Zhang   } else {
1201a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1202a4e85ca8SHong Zhang   }
1203a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1204a4e85ca8SHong Zhang }
1205a4e85ca8SHong Zhang 
120624121865SAdrian Maldonado 
120724121865SAdrian Maldonado /* Creates a GlobalToLocal mapping with a Local and Global section. This is akin to the routine DMGetLocalToGlobalMapping but without the need of providing a dm.
120824121865SAdrian Maldonado */
120924121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
121024121865SAdrian Maldonado {
121124121865SAdrian Maldonado   PetscErrorCode ierr;
121224121865SAdrian Maldonado   PetscInt       i, size, dof;
121324121865SAdrian Maldonado   PetscInt       *glob2loc;
121424121865SAdrian Maldonado 
121524121865SAdrian Maldonado   PetscFunctionBegin;
121624121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
121724121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
121824121865SAdrian Maldonado 
121924121865SAdrian Maldonado   for (i = 0; i < size; i++) {
122024121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
122124121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
122224121865SAdrian Maldonado     glob2loc[i] = dof;
122324121865SAdrian Maldonado   }
122424121865SAdrian Maldonado 
122524121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
122624121865SAdrian Maldonado #if 0
122724121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
122824121865SAdrian Maldonado #endif
122924121865SAdrian Maldonado   PetscFunctionReturn(0);
123024121865SAdrian Maldonado }
123124121865SAdrian Maldonado 
12321ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
12331ad426b7SHong Zhang {
12341ad426b7SHong Zhang   PetscErrorCode ierr;
123524121865SAdrian Maldonado   PetscMPIInt    rank, size;
12361ad426b7SHong Zhang   DM_Network     *network = (DM_Network*) dm->data;
1237a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1238840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
123924121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1240a4e85ca8SHong Zhang   Mat            Juser;
1241bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1242447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1243a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
12445cf7da58SHong Zhang   MPI_Comm       comm;
124524121865SAdrian Maldonado   MatType        mtype;
12465cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
12475cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
12485cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
12491ad426b7SHong Zhang 
12501ad426b7SHong Zhang   PetscFunctionBegin;
125124121865SAdrian Maldonado   mtype = dm->mattype;
125224121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
125324121865SAdrian Maldonado 
125424121865SAdrian Maldonado   if (isNest) {
12550731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
125624121865SAdrian Maldonado     PetscInt   eDof, vDof;
125724121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
125824121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
125924121865SAdrian Maldonado 
126024121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
126124121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
126224121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
126324121865SAdrian Maldonado 
126424121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
126524121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
126624121865SAdrian Maldonado 
126724121865SAdrian Maldonado     ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr);
126824121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
126924121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
127024121865SAdrian Maldonado 
127124121865SAdrian Maldonado     ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr);
127224121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
127324121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
127424121865SAdrian Maldonado 
127524121865SAdrian Maldonado     ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr);
127624121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
127724121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
127824121865SAdrian Maldonado 
127924121865SAdrian Maldonado     ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr);
128024121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
128124121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
128224121865SAdrian Maldonado 
128324121865SAdrian Maldonado     bA[0][0] = j11;
128424121865SAdrian Maldonado     bA[0][1] = j12;
128524121865SAdrian Maldonado     bA[1][0] = j21;
128624121865SAdrian Maldonado     bA[1][1] = j22;
128724121865SAdrian Maldonado 
128824121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
128924121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
129024121865SAdrian Maldonado 
129124121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
129224121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
129324121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
129424121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
129524121865SAdrian Maldonado 
129624121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
129724121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
129824121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
129924121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
130024121865SAdrian Maldonado 
130124121865SAdrian Maldonado     ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
130224121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
130324121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
130424121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
130524121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
130624121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
130724121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
130824121865SAdrian Maldonado 
130924121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
131024121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
131124121865SAdrian Maldonado 
131224121865SAdrian Maldonado     /* Free structures */
131324121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
131424121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
131524121865SAdrian Maldonado 
131624121865SAdrian Maldonado     PetscFunctionReturn(0);
131724121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1318a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1319bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1320bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
13211ad426b7SHong Zhang     PetscFunctionReturn(0);
13221ad426b7SHong Zhang   }
13231ad426b7SHong Zhang 
1324bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
13252a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1326bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1327bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
13282a945128SHong Zhang 
13292a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
13302a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
133189898e50SHong Zhang 
133289898e50SHong Zhang   /* (1) Set matrix preallocation */
133389898e50SHong Zhang   /*------------------------------*/
1334840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1335840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1336840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1337840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1338840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1339840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1340840c2264SHong Zhang 
134189898e50SHong Zhang   /* Set preallocation for edges */
134289898e50SHong Zhang   /*-----------------------------*/
1343840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1344840c2264SHong Zhang 
1345bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1346840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1347840c2264SHong Zhang     /* Get row indices */
1348840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1349840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1350840c2264SHong Zhang     if (nrows) {
1351840c2264SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1352840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1353840c2264SHong Zhang 
13545cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1355840c2264SHong Zhang       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1356840c2264SHong Zhang       for (v=0; v<2; v++) {
1357840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1358840c2264SHong Zhang 
1359840c2264SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1360840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
13615cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1362840c2264SHong Zhang       }
1363840c2264SHong Zhang 
136489898e50SHong Zhang       /* Set preallocation for edge self */
1365840c2264SHong Zhang       cstart = rstart;
1366840c2264SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
13675cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1368840c2264SHong Zhang     }
1369840c2264SHong Zhang   }
1370840c2264SHong Zhang 
137189898e50SHong Zhang   /* Set preallocation for vertices */
137289898e50SHong Zhang   /*--------------------------------*/
1373840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1374840c2264SHong Zhang   if (vEnd - vStart) {
1375840c2264SHong Zhang     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1376840c2264SHong Zhang     vptr = network->Jvptr;
1377840c2264SHong Zhang     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1378840c2264SHong Zhang   }
1379840c2264SHong Zhang 
1380840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1381840c2264SHong Zhang     /* Get row indices */
1382840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1383840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1384840c2264SHong Zhang     if (!nrows) continue;
1385840c2264SHong Zhang 
1386bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1387bdcb62a2SHong Zhang     if (ghost) {
1388bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1389bdcb62a2SHong Zhang     } else {
1390bdcb62a2SHong Zhang       rows_v = rows;
1391bdcb62a2SHong Zhang     }
1392bdcb62a2SHong Zhang 
1393bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1394840c2264SHong Zhang 
1395840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1396840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1397840c2264SHong Zhang 
1398840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1399840c2264SHong Zhang       /* Supporting edges */
1400840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1401840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1402840c2264SHong Zhang 
1403840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1404bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1405840c2264SHong Zhang 
1406840c2264SHong Zhang       /* Connected vertices */
1407840c2264SHong Zhang       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
1408840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1409840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1410840c2264SHong Zhang 
1411840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1412840c2264SHong Zhang 
1413840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1414e102a522SHong Zhang       if (ghost_vc||ghost) {
1415e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1416e102a522SHong Zhang       } else {
1417e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1418e102a522SHong Zhang       }
1419e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1420840c2264SHong Zhang     }
1421840c2264SHong Zhang 
142289898e50SHong Zhang     /* Set preallocation for vertex self */
1423840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1424840c2264SHong Zhang     if (!ghost) {
1425840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1426840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1427bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1428840c2264SHong Zhang     }
1429bdcb62a2SHong Zhang     if (ghost) {
1430bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1431bdcb62a2SHong Zhang     }
1432840c2264SHong Zhang   }
1433840c2264SHong Zhang 
1434840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1435840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
14365cf7da58SHong Zhang 
14375cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
14385cf7da58SHong Zhang 
14395cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1440840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1441840c2264SHong Zhang 
1442840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1443840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1444840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1445e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1446e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1447840c2264SHong Zhang   }
1448840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1449840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1450840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1451840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1452840c2264SHong Zhang 
14535cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
14545cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
14555cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
14565cf7da58SHong Zhang 
14575cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
14585cf7da58SHong Zhang 
145989898e50SHong Zhang   /* (2) Set matrix entries for edges */
146089898e50SHong Zhang   /*----------------------------------*/
14611ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1462bfbc38dcSHong Zhang     /* Get row indices */
14631ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
146417df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
14654b976069SHong Zhang     if (nrows) {
14664b976069SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1467bdcb62a2SHong Zhang 
146817df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
14691ad426b7SHong Zhang 
1470bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
147117df6e9eSHong Zhang       ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr);
1472bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1473bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1474883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
14753e97b6e8SHong Zhang 
1476a4e85ca8SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1477a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1478bfbc38dcSHong Zhang       }
147917df6e9eSHong Zhang 
1480bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
14813e97b6e8SHong Zhang       cstart = rstart;
1482a4e85ca8SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1483a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
14841ad426b7SHong Zhang     }
14854b976069SHong Zhang   }
14861ad426b7SHong Zhang 
1487bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
148883b2e829SHong Zhang   /*---------------------------------*/
14891ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1490bfbc38dcSHong Zhang     /* Get row indices */
1491596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1492596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
14934b976069SHong Zhang     if (!nrows) continue;
1494596e729fSHong Zhang 
1495bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1496bdcb62a2SHong Zhang     if (ghost) {
1497bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1498bdcb62a2SHong Zhang     } else {
1499bdcb62a2SHong Zhang       rows_v = rows;
1500bdcb62a2SHong Zhang     }
1501bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1502596e729fSHong Zhang 
1503bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1504596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1505596e729fSHong Zhang 
1506596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1507bfbc38dcSHong Zhang       /* Supporting edges */
1508596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1509596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1510596e729fSHong Zhang 
1511a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1512bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1513596e729fSHong Zhang 
1514bfbc38dcSHong Zhang       /* Connected vertices */
151544aca652SHong Zhang       ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr);
15162a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
15172a945128SHong Zhang 
151844aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
151944aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1520a4e85ca8SHong Zhang 
1521a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1522bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1523596e729fSHong Zhang     }
1524596e729fSHong Zhang 
1525bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
15261ad426b7SHong Zhang     if (!ghost) {
1527596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1528a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1529bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1530bdcb62a2SHong Zhang     }
1531bdcb62a2SHong Zhang     if (ghost) {
1532bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1533bdcb62a2SHong Zhang     }
15341ad426b7SHong Zhang   }
1535a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1536bdcb62a2SHong Zhang 
15371ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15381ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1539dd6f46cdSHong Zhang 
15405f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
15415f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
15425f2c45f1SShri Abhyankar }
15435f2c45f1SShri Abhyankar 
15445f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
15455f2c45f1SShri Abhyankar {
15465f2c45f1SShri Abhyankar   PetscErrorCode ierr;
15475f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
15485f2c45f1SShri Abhyankar 
15495f2c45f1SShri Abhyankar   PetscFunctionBegin;
15508415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
155183b2e829SHong Zhang   if (network->Je) {
155283b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
155383b2e829SHong Zhang   }
155483b2e829SHong Zhang   if (network->Jv) {
1555883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
155683b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
15571ad426b7SHong Zhang   }
155813c2a604SAdrian Maldonado 
155913c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
156013c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
156113c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
156213c2a604SAdrian Maldonado   if (network->vertex.sf) {
156313c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
156413c2a604SAdrian Maldonado   }
156513c2a604SAdrian Maldonado   /* edge */
156613c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
156713c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
156813c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
156913c2a604SAdrian Maldonado   if (network->edge.sf) {
157013c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
157113c2a604SAdrian Maldonado   }
15725f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
15735f2c45f1SShri Abhyankar   network->edges = NULL;
15745f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
15755f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
157683b2e829SHong Zhang 
15775f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
15785f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
15795f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
15805f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
15815f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
15825f2c45f1SShri Abhyankar }
15835f2c45f1SShri Abhyankar 
15845f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
15855f2c45f1SShri Abhyankar {
15865f2c45f1SShri Abhyankar   PetscErrorCode ierr;
15875f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
15885f2c45f1SShri Abhyankar 
15895f2c45f1SShri Abhyankar   PetscFunctionBegin;
15905f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
15915f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
15925f2c45f1SShri Abhyankar }
15935f2c45f1SShri Abhyankar 
15945f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
15955f2c45f1SShri Abhyankar {
15965f2c45f1SShri Abhyankar   PetscErrorCode ierr;
15975f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
15985f2c45f1SShri Abhyankar 
15995f2c45f1SShri Abhyankar   PetscFunctionBegin;
16005f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
16015f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16025f2c45f1SShri Abhyankar }
16035f2c45f1SShri Abhyankar 
16045f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
16055f2c45f1SShri Abhyankar {
16065f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16075f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16085f2c45f1SShri Abhyankar 
16095f2c45f1SShri Abhyankar   PetscFunctionBegin;
16105f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
16115f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16125f2c45f1SShri Abhyankar }
16135f2c45f1SShri Abhyankar 
16145f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
16155f2c45f1SShri Abhyankar {
16165f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16175f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16185f2c45f1SShri Abhyankar 
16195f2c45f1SShri Abhyankar   PetscFunctionBegin;
16205f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
16215f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16225f2c45f1SShri Abhyankar }
16235f2c45f1SShri Abhyankar 
16245f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
16255f2c45f1SShri Abhyankar {
16265f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16275f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16285f2c45f1SShri Abhyankar 
16295f2c45f1SShri Abhyankar   PetscFunctionBegin;
16305f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
16315f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16325f2c45f1SShri Abhyankar }
1633