xref: /petsc/src/dm/impls/network/network.c (revision 6fefedf4bbd7aa4e9e6bb83314d11516b7b0319e)
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);
60*6fefedf4SHong Zhang   if ((network->nVertices >= 0 || network->NVertices >= 0) && (network->nVertices != nV || network->NVertices != 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->nVertices,network->NVertices);
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   }
67*6fefedf4SHong Zhang   network->nVertices = nV;
68*6fefedf4SHong Zhang   network->NVertices = 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;
129*6fefedf4SHong Zhang   if (network->nVertices) {
130*6fefedf4SHong Zhang     ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr);
1315f2c45f1SShri Abhyankar   }
132*6fefedf4SHong Zhang   ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr);
133*6fefedf4SHong Zhang   if (network->nVertices) {
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++) {
1507b6afd5bSHong Zhang     if (i < network->vStart) {
151e85e6aecSHong Zhang       network->header[i].index = i;
1527b6afd5bSHong Zhang     } else {
153e85e6aecSHong Zhang       network->header[i].index = i - network->vStart;
1547b6afd5bSHong 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 
2617b6afd5bSHong Zhang /*@
262e85e6aecSHong Zhang   DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge.
2637b6afd5bSHong Zhang 
2647b6afd5bSHong Zhang   Not Collective
2657b6afd5bSHong Zhang 
2667b6afd5bSHong Zhang   Input Parameters:
2677b6afd5bSHong Zhang + dm - DMNetwork object
268e85e6aecSHong Zhang - p  - edge point
2697b6afd5bSHong Zhang 
2707b6afd5bSHong Zhang   Output Paramters:
271e85e6aecSHong Zhang . index - user global numbering for the edge
2727b6afd5bSHong Zhang 
2737b6afd5bSHong Zhang   Level: intermediate
2747b6afd5bSHong Zhang 
275e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex
2767b6afd5bSHong Zhang @*/
277e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index)
2787b6afd5bSHong Zhang {
2797b6afd5bSHong Zhang   PetscErrorCode    ierr;
2807b6afd5bSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
2817b6afd5bSHong Zhang   PetscInt          offsetp;
2827b6afd5bSHong Zhang   DMNetworkComponentHeader header;
2837b6afd5bSHong Zhang 
2847b6afd5bSHong Zhang   PetscFunctionBegin;
2857b6afd5bSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
2867b6afd5bSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
287e85e6aecSHong Zhang   *index = header->index;
2887b6afd5bSHong Zhang   PetscFunctionReturn(0);
2897b6afd5bSHong Zhang }
2907b6afd5bSHong Zhang 
2915f2c45f1SShri Abhyankar /*@
292e85e6aecSHong Zhang   DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex.
293e85e6aecSHong Zhang 
294e85e6aecSHong Zhang   Not Collective
295e85e6aecSHong Zhang 
296e85e6aecSHong Zhang   Input Parameters:
297e85e6aecSHong Zhang + dm - DMNetwork object
298e85e6aecSHong Zhang - p  - vertex point
299e85e6aecSHong Zhang 
300e85e6aecSHong Zhang   Output Paramters:
301e85e6aecSHong Zhang . index - user global numbering for the vertex
302e85e6aecSHong Zhang 
303e85e6aecSHong Zhang   Level: intermediate
304e85e6aecSHong Zhang 
305e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex
306e85e6aecSHong Zhang @*/
307e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index)
308e85e6aecSHong Zhang {
309e85e6aecSHong Zhang   PetscErrorCode    ierr;
310e85e6aecSHong Zhang   DM_Network        *network = (DM_Network*)dm->data;
311e85e6aecSHong Zhang   PetscInt          offsetp;
312e85e6aecSHong Zhang   DMNetworkComponentHeader header;
313e85e6aecSHong Zhang 
314e85e6aecSHong Zhang   PetscFunctionBegin;
315e85e6aecSHong Zhang   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
316e85e6aecSHong Zhang   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
317e85e6aecSHong Zhang   *index = header->index;
318e85e6aecSHong Zhang   PetscFunctionReturn(0);
319e85e6aecSHong Zhang }
320e85e6aecSHong Zhang 
321e85e6aecSHong Zhang /*@
322325661f6SSatish Balay   DMNetworkAddComponent - Adds a network component at the given point (vertex/edge)
3235f2c45f1SShri Abhyankar 
3245f2c45f1SShri Abhyankar   Not Collective
3255f2c45f1SShri Abhyankar 
3265f2c45f1SShri Abhyankar   Input Parameters:
3275f2c45f1SShri Abhyankar + dm           - The DMNetwork object
3285f2c45f1SShri Abhyankar . p            - vertex/edge point
3295f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component
3305f2c45f1SShri Abhyankar - compvalue    - pointer to the data structure for the component
3315f2c45f1SShri Abhyankar 
3325f2c45f1SShri Abhyankar   Level: intermediate
3335f2c45f1SShri Abhyankar 
3345f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent
3355f2c45f1SShri Abhyankar @*/
3365f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue)
3375f2c45f1SShri Abhyankar {
3385f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
33943a39a44SBarry Smith   DMNetworkComponent       *component = &network->component[componentkey];
3405f2c45f1SShri Abhyankar   DMNetworkComponentHeader header = &network->header[p];
3415f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue = &network->cvalue[p];
3425f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
3435f2c45f1SShri Abhyankar 
3445f2c45f1SShri Abhyankar   PetscFunctionBegin;
345fa58f0a9SHong 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);
346fa58f0a9SHong Zhang 
34743a39a44SBarry Smith   header->size[header->ndata] = component->size;
34843a39a44SBarry Smith   ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr);
3495f2c45f1SShri Abhyankar   header->key[header->ndata] = componentkey;
3505f2c45f1SShri Abhyankar   if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1];
3515f2c45f1SShri Abhyankar 
3525f2c45f1SShri Abhyankar   cvalue->data[header->ndata] = (void*)compvalue;
3535f2c45f1SShri Abhyankar   header->ndata++;
3545f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3555f2c45f1SShri Abhyankar }
3565f2c45f1SShri Abhyankar 
3575f2c45f1SShri Abhyankar /*@
3585f2c45f1SShri Abhyankar   DMNetworkGetNumComponents - Get the number of components at a vertex/edge
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 
3665f2c45f1SShri Abhyankar   Output Parameters:
3675f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge
3685f2c45f1SShri Abhyankar 
3695f2c45f1SShri Abhyankar   Level: intermediate
3705f2c45f1SShri Abhyankar 
3715f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent
3725f2c45f1SShri Abhyankar @*/
3735f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents)
3745f2c45f1SShri Abhyankar {
3755f2c45f1SShri Abhyankar   PetscErrorCode ierr;
3765f2c45f1SShri Abhyankar   PetscInt       offset;
3775f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
3785f2c45f1SShri Abhyankar 
3795f2c45f1SShri Abhyankar   PetscFunctionBegin;
3805f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr);
3815f2c45f1SShri Abhyankar   *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata;
3825f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
3835f2c45f1SShri Abhyankar }
3845f2c45f1SShri Abhyankar 
3855f2c45f1SShri Abhyankar /*@
386a730d845SHong Zhang   DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the
3875f2c45f1SShri Abhyankar                                     component value from the component data array
3885f2c45f1SShri Abhyankar 
3895f2c45f1SShri Abhyankar   Not Collective
3905f2c45f1SShri Abhyankar 
3915f2c45f1SShri Abhyankar   Input Parameters:
3925f2c45f1SShri Abhyankar + dm      - The DMNetwork object
3935f2c45f1SShri Abhyankar . p       - vertex/edge point
3945f2c45f1SShri Abhyankar - compnum - component number
3955f2c45f1SShri Abhyankar 
3965f2c45f1SShri Abhyankar   Output Parameters:
3975f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component
3985f2c45f1SShri Abhyankar - offset  - offset into the component data array associated with the vertex/edge point
3995f2c45f1SShri Abhyankar 
4005f2c45f1SShri Abhyankar   Notes:
4015f2c45f1SShri Abhyankar   Typical usage:
4025f2c45f1SShri Abhyankar 
4035f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray(dm, &arr);
4045f2c45f1SShri Abhyankar   DMNetworkGetVertex/EdgeRange(dm,&Start,&End);
4055f2c45f1SShri Abhyankar   Loop over vertices or edges
4065f2c45f1SShri Abhyankar     DMNetworkGetNumComponents(dm,v,&numcomps);
4075f2c45f1SShri Abhyankar     Loop over numcomps
408a730d845SHong Zhang       DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset);
4095f2c45f1SShri Abhyankar       compdata = (UserCompDataType)(arr+offset);
4105f2c45f1SShri Abhyankar 
4115f2c45f1SShri Abhyankar   Level: intermediate
4125f2c45f1SShri Abhyankar 
4135f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray,
4145f2c45f1SShri Abhyankar @*/
415a730d845SHong Zhang PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset)
4165f2c45f1SShri Abhyankar {
4175f2c45f1SShri Abhyankar   PetscErrorCode           ierr;
4185f2c45f1SShri Abhyankar   PetscInt                 offsetp;
4195f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
4205f2c45f1SShri Abhyankar   DM_Network               *network = (DM_Network*)dm->data;
4215f2c45f1SShri Abhyankar 
4225f2c45f1SShri Abhyankar   PetscFunctionBegin;
4235f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
4245f2c45f1SShri Abhyankar   header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp);
425d36f4e81SHong Zhang   if (compkey) *compkey = header->key[compnum];
426d36f4e81SHong Zhang   if (offset) *offset  = offsetp+network->dataheadersize+header->offset[compnum];
4275f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4285f2c45f1SShri Abhyankar }
4295f2c45f1SShri Abhyankar 
4305f2c45f1SShri Abhyankar /*@
4315f2c45f1SShri Abhyankar   DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector.
4325f2c45f1SShri Abhyankar 
4335f2c45f1SShri Abhyankar   Not Collective
4345f2c45f1SShri Abhyankar 
4355f2c45f1SShri Abhyankar   Input Parameters:
4365f2c45f1SShri Abhyankar + dm     - The DMNetwork object
4375f2c45f1SShri Abhyankar - p      - the edge/vertex point
4385f2c45f1SShri Abhyankar 
4395f2c45f1SShri Abhyankar   Output Parameters:
4405f2c45f1SShri Abhyankar . offset - the offset
4415f2c45f1SShri Abhyankar 
4425f2c45f1SShri Abhyankar   Level: intermediate
4435f2c45f1SShri Abhyankar 
4445f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
4455f2c45f1SShri Abhyankar @*/
4465f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset)
4475f2c45f1SShri Abhyankar {
4485f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4495f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4505f2c45f1SShri Abhyankar 
4515f2c45f1SShri Abhyankar   PetscFunctionBegin;
4525f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr);
4535f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4545f2c45f1SShri Abhyankar }
4555f2c45f1SShri Abhyankar 
4565f2c45f1SShri Abhyankar /*@
4575f2c45f1SShri Abhyankar   DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector.
4585f2c45f1SShri Abhyankar 
4595f2c45f1SShri Abhyankar   Not Collective
4605f2c45f1SShri Abhyankar 
4615f2c45f1SShri Abhyankar   Input Parameters:
4625f2c45f1SShri Abhyankar + dm      - The DMNetwork object
4635f2c45f1SShri Abhyankar - p       - the edge/vertex point
4645f2c45f1SShri Abhyankar 
4655f2c45f1SShri Abhyankar   Output Parameters:
4665f2c45f1SShri Abhyankar . offsetg - the offset
4675f2c45f1SShri Abhyankar 
4685f2c45f1SShri Abhyankar   Level: intermediate
4695f2c45f1SShri Abhyankar 
4705f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector
4715f2c45f1SShri Abhyankar @*/
4725f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg)
4735f2c45f1SShri Abhyankar {
4745f2c45f1SShri Abhyankar   PetscErrorCode ierr;
4755f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
4765f2c45f1SShri Abhyankar 
4775f2c45f1SShri Abhyankar   PetscFunctionBegin;
4785f78ed8bSShri Abhyankar   ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr);
479*6fefedf4SHong Zhang   if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */
4805f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
4815f2c45f1SShri Abhyankar }
4825f2c45f1SShri Abhyankar 
48324121865SAdrian Maldonado /*@
48424121865SAdrian Maldonado   DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector.
48524121865SAdrian Maldonado 
48624121865SAdrian Maldonado   Not Collective
48724121865SAdrian Maldonado 
48824121865SAdrian Maldonado   Input Parameters:
48924121865SAdrian Maldonado + dm     - The DMNetwork object
49024121865SAdrian Maldonado - p      - the edge point
49124121865SAdrian Maldonado 
49224121865SAdrian Maldonado   Output Parameters:
49324121865SAdrian Maldonado . offset - the offset
49424121865SAdrian Maldonado 
49524121865SAdrian Maldonado   Level: intermediate
49624121865SAdrian Maldonado 
49724121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
49824121865SAdrian Maldonado @*/
49924121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset)
50024121865SAdrian Maldonado {
50124121865SAdrian Maldonado   PetscErrorCode ierr;
50224121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
50324121865SAdrian Maldonado 
50424121865SAdrian Maldonado   PetscFunctionBegin;
50524121865SAdrian Maldonado 
50624121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr);
50724121865SAdrian Maldonado   PetscFunctionReturn(0);
50824121865SAdrian Maldonado }
50924121865SAdrian Maldonado 
51024121865SAdrian Maldonado /*@
51124121865SAdrian Maldonado   DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector.
51224121865SAdrian Maldonado 
51324121865SAdrian Maldonado   Not Collective
51424121865SAdrian Maldonado 
51524121865SAdrian Maldonado   Input Parameters:
51624121865SAdrian Maldonado + dm     - The DMNetwork object
51724121865SAdrian Maldonado - p      - the vertex point
51824121865SAdrian Maldonado 
51924121865SAdrian Maldonado   Output Parameters:
52024121865SAdrian Maldonado . offset - the offset
52124121865SAdrian Maldonado 
52224121865SAdrian Maldonado   Level: intermediate
52324121865SAdrian Maldonado 
52424121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector
52524121865SAdrian Maldonado @*/
52624121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset)
52724121865SAdrian Maldonado {
52824121865SAdrian Maldonado   PetscErrorCode ierr;
52924121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
53024121865SAdrian Maldonado 
53124121865SAdrian Maldonado   PetscFunctionBegin;
53224121865SAdrian Maldonado 
53324121865SAdrian Maldonado   p -= network->vStart;
53424121865SAdrian Maldonado 
53524121865SAdrian Maldonado   ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr);
53624121865SAdrian Maldonado   PetscFunctionReturn(0);
53724121865SAdrian Maldonado }
5385f2c45f1SShri Abhyankar /*@
5395f2c45f1SShri Abhyankar   DMNetworkAddNumVariables - Add number of variables associated with a given point.
5405f2c45f1SShri Abhyankar 
5415f2c45f1SShri Abhyankar   Not Collective
5425f2c45f1SShri Abhyankar 
5435f2c45f1SShri Abhyankar   Input Parameters:
5445f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
5455f2c45f1SShri Abhyankar . p    - the vertex/edge point
5465f2c45f1SShri Abhyankar - nvar - number of additional variables
5475f2c45f1SShri Abhyankar 
5485f2c45f1SShri Abhyankar   Level: intermediate
5495f2c45f1SShri Abhyankar 
5505f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables
5515f2c45f1SShri Abhyankar @*/
5525f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar)
5535f2c45f1SShri Abhyankar {
5545f2c45f1SShri Abhyankar   PetscErrorCode ierr;
5555f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
5565f2c45f1SShri Abhyankar 
5575f2c45f1SShri Abhyankar   PetscFunctionBegin;
5585f2c45f1SShri Abhyankar   ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr);
5595f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
5605f2c45f1SShri Abhyankar }
5615f2c45f1SShri Abhyankar 
56227f51fceSHong Zhang /*@
56327f51fceSHong Zhang   DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point.
56427f51fceSHong Zhang 
56527f51fceSHong Zhang   Not Collective
56627f51fceSHong Zhang 
56727f51fceSHong Zhang   Input Parameters:
56827f51fceSHong Zhang + dm   - The DMNetworkObject
56927f51fceSHong Zhang - p    - the vertex/edge point
57027f51fceSHong Zhang 
57127f51fceSHong Zhang   Output Parameters:
57227f51fceSHong Zhang . nvar - number of variables
57327f51fceSHong Zhang 
57427f51fceSHong Zhang   Level: intermediate
57527f51fceSHong Zhang 
57627f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables
57727f51fceSHong Zhang @*/
57827f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar)
57927f51fceSHong Zhang {
58027f51fceSHong Zhang   PetscErrorCode ierr;
58127f51fceSHong Zhang   DM_Network     *network = (DM_Network*)dm->data;
58227f51fceSHong Zhang 
58327f51fceSHong Zhang   PetscFunctionBegin;
58427f51fceSHong Zhang   ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
58527f51fceSHong Zhang   PetscFunctionReturn(0);
58627f51fceSHong Zhang }
58727f51fceSHong Zhang 
5885f2c45f1SShri Abhyankar /*@
5895f2c45f1SShri Abhyankar   DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point.
5905f2c45f1SShri Abhyankar 
5915f2c45f1SShri Abhyankar   Not Collective
5925f2c45f1SShri Abhyankar 
5935f2c45f1SShri Abhyankar   Input Parameters:
5945f2c45f1SShri Abhyankar + dm   - The DMNetworkObject
5955f2c45f1SShri Abhyankar . p    - the vertex/edge point
5965f2c45f1SShri Abhyankar - nvar - number of variables
5975f2c45f1SShri Abhyankar 
5985f2c45f1SShri Abhyankar   Level: intermediate
5995f2c45f1SShri Abhyankar 
6005f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables
6015f2c45f1SShri Abhyankar @*/
6025f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar)
6035f2c45f1SShri Abhyankar {
6045f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6055f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6065f2c45f1SShri Abhyankar 
6075f2c45f1SShri Abhyankar   PetscFunctionBegin;
6085f2c45f1SShri Abhyankar   ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr);
6095f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6105f2c45f1SShri Abhyankar }
6115f2c45f1SShri Abhyankar 
6125f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This
6135f2c45f1SShri Abhyankar    function is called during DMSetUp() */
6145f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm)
6155f2c45f1SShri Abhyankar {
6165f2c45f1SShri Abhyankar   PetscErrorCode              ierr;
6175f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6185f2c45f1SShri Abhyankar   PetscInt                    arr_size;
6195f2c45f1SShri Abhyankar   PetscInt                    p,offset,offsetp;
6205f2c45f1SShri Abhyankar   DMNetworkComponentHeader header;
6215f2c45f1SShri Abhyankar   DMNetworkComponentValue  cvalue;
6225f2c45f1SShri Abhyankar   DMNetworkComponentGenericDataType      *componentdataarray;
6235f2c45f1SShri Abhyankar   PetscInt ncomp, i;
6245f2c45f1SShri Abhyankar 
6255f2c45f1SShri Abhyankar   PetscFunctionBegin;
6265f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr);
6275f2c45f1SShri Abhyankar   ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr);
62875b160a0SShri Abhyankar   ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr);
6295f2c45f1SShri Abhyankar   componentdataarray = network->componentdataarray;
6305f2c45f1SShri Abhyankar   for (p = network->pStart; p < network->pEnd; p++) {
6315f2c45f1SShri Abhyankar     ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr);
6325f2c45f1SShri Abhyankar     /* Copy header */
6335f2c45f1SShri Abhyankar     header = &network->header[p];
634302440fdSBarry Smith     ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
6355f2c45f1SShri Abhyankar     /* Copy data */
6365f2c45f1SShri Abhyankar     cvalue = &network->cvalue[p];
6375f2c45f1SShri Abhyankar     ncomp = header->ndata;
6385f2c45f1SShri Abhyankar     for (i = 0; i < ncomp; i++) {
6395f2c45f1SShri Abhyankar       offset = offsetp + network->dataheadersize + header->offset[i];
640302440fdSBarry Smith       ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr);
6415f2c45f1SShri Abhyankar     }
6425f2c45f1SShri Abhyankar   }
6435f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6445f2c45f1SShri Abhyankar }
6455f2c45f1SShri Abhyankar 
6465f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */
6475f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm)
6485f2c45f1SShri Abhyankar {
6495f2c45f1SShri Abhyankar   PetscErrorCode ierr;
6505f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6515f2c45f1SShri Abhyankar 
6525f2c45f1SShri Abhyankar   PetscFunctionBegin;
6535f2c45f1SShri Abhyankar   ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr);
6545f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6555f2c45f1SShri Abhyankar }
6565f2c45f1SShri Abhyankar 
6575f2c45f1SShri Abhyankar /*@C
6585f2c45f1SShri Abhyankar   DMNetworkGetComponentDataArray - Returns the component data array
6595f2c45f1SShri Abhyankar 
6605f2c45f1SShri Abhyankar   Not Collective
6615f2c45f1SShri Abhyankar 
6625f2c45f1SShri Abhyankar   Input Parameters:
6635f2c45f1SShri Abhyankar . dm - The DMNetwork Object
6645f2c45f1SShri Abhyankar 
6655f2c45f1SShri Abhyankar   Output Parameters:
6665f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components
6675f2c45f1SShri Abhyankar 
6685f2c45f1SShri Abhyankar   Level: intermediate
6695f2c45f1SShri Abhyankar 
670a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents
6715f2c45f1SShri Abhyankar @*/
6725f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray)
6735f2c45f1SShri Abhyankar {
6745f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
6755f2c45f1SShri Abhyankar 
6765f2c45f1SShri Abhyankar   PetscFunctionBegin;
6775f2c45f1SShri Abhyankar   *componentdataarray = network->componentdataarray;
6785f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
6795f2c45f1SShri Abhyankar }
6805f2c45f1SShri Abhyankar 
68124121865SAdrian Maldonado /* Get a subsection from a range of points */
68224121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection)
68324121865SAdrian Maldonado {
68424121865SAdrian Maldonado   PetscErrorCode ierr;
68524121865SAdrian Maldonado   PetscInt       i, nvar;
68624121865SAdrian Maldonado 
68724121865SAdrian Maldonado   PetscFunctionBegin;
68824121865SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr);
68924121865SAdrian Maldonado   ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr);
69024121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
69124121865SAdrian Maldonado     ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr);
69224121865SAdrian Maldonado     ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr);
69324121865SAdrian Maldonado   }
69424121865SAdrian Maldonado 
69524121865SAdrian Maldonado   ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr);
69624121865SAdrian Maldonado   PetscFunctionReturn(0);
69724121865SAdrian Maldonado }
69824121865SAdrian Maldonado 
69924121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */
70024121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map)
70124121865SAdrian Maldonado {
70224121865SAdrian Maldonado   PetscErrorCode ierr;
70324121865SAdrian Maldonado   PetscInt       i, *subpoints;
70424121865SAdrian Maldonado 
70524121865SAdrian Maldonado   PetscFunctionBegin;
70624121865SAdrian Maldonado   /* Create index sets to map from "points" to "subpoints" */
70724121865SAdrian Maldonado   ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr);
70824121865SAdrian Maldonado   for (i = pstart; i < pend; i++) {
70924121865SAdrian Maldonado     subpoints[i - pstart] = i;
71024121865SAdrian Maldonado   }
711459726d8SSatish Balay   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr);
71224121865SAdrian Maldonado   ierr = PetscFree(subpoints);CHKERRQ(ierr);
71324121865SAdrian Maldonado   PetscFunctionReturn(0);
71424121865SAdrian Maldonado }
71524121865SAdrian Maldonado 
71624121865SAdrian Maldonado /*@
71724121865SAdrian Maldonado   DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute.
71824121865SAdrian Maldonado 
71924121865SAdrian Maldonado   Collective
72024121865SAdrian Maldonado 
72124121865SAdrian Maldonado   Input Parameters:
72224121865SAdrian Maldonado . dm   - The DMNetworkObject
72324121865SAdrian Maldonado 
72424121865SAdrian Maldonado   Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are:
72524121865SAdrian Maldonado 
72624121865SAdrian Maldonado   points = [0 1 2 3 4 5 6]
72724121865SAdrian Maldonado 
72824121865SAdrian Maldonado   where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]).
72924121865SAdrian Maldonado 
73024121865SAdrian Maldonado   With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset.
73124121865SAdrian Maldonado 
73224121865SAdrian Maldonado   Level: intermediate
73324121865SAdrian Maldonado 
73424121865SAdrian Maldonado @*/
73524121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm)
73624121865SAdrian Maldonado {
73724121865SAdrian Maldonado   PetscErrorCode ierr;
73824121865SAdrian Maldonado   MPI_Comm       comm;
7399852e123SBarry Smith   PetscMPIInt    rank, size;
74024121865SAdrian Maldonado   DM_Network     *network = (DM_Network*)dm->data;
74124121865SAdrian Maldonado 
742eab1376dSHong Zhang   PetscFunctionBegin;
74324121865SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
74424121865SAdrian Maldonado   ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr);
7459852e123SBarry Smith   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
74624121865SAdrian Maldonado 
74724121865SAdrian Maldonado   /* Create maps for vertices and edges */
74824121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr);
74924121865SAdrian Maldonado   ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr);
75024121865SAdrian Maldonado 
75124121865SAdrian Maldonado   /* Create local sub-sections */
75224121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr);
75324121865SAdrian Maldonado   ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr);
75424121865SAdrian Maldonado 
7559852e123SBarry Smith   if (size > 1) {
75624121865SAdrian Maldonado     ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr);
75724121865SAdrian Maldonado     ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr);
75824121865SAdrian Maldonado   ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr);
75924121865SAdrian Maldonado   ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr);
76024121865SAdrian Maldonado   } else {
76124121865SAdrian Maldonado   /* create structures for vertex */
76224121865SAdrian Maldonado   ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr);
76324121865SAdrian Maldonado   /* create structures for edge */
76424121865SAdrian Maldonado   ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr);
76524121865SAdrian Maldonado   }
76624121865SAdrian Maldonado 
76724121865SAdrian Maldonado 
76824121865SAdrian Maldonado   /* Add viewers */
76924121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr);
77024121865SAdrian Maldonado   ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr);
77124121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr);
77224121865SAdrian Maldonado   ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr);
77324121865SAdrian Maldonado 
77424121865SAdrian Maldonado   PetscFunctionReturn(0);
77524121865SAdrian Maldonado }
7767b6afd5bSHong Zhang 
7775f2c45f1SShri Abhyankar /*@
7785f2c45f1SShri Abhyankar   DMNetworkDistribute - Distributes the network and moves associated component data.
7795f2c45f1SShri Abhyankar 
7805f2c45f1SShri Abhyankar   Collective
7815f2c45f1SShri Abhyankar 
7825f2c45f1SShri Abhyankar   Input Parameter:
783d3464fd4SAdrian Maldonado + DM - the DMNetwork object
7845f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default
7855f2c45f1SShri Abhyankar 
7865f2c45f1SShri Abhyankar   Notes:
7878b171c8eSHong Zhang   Distributes the network with <overlap>-overlapping partitioning of the edges.
7885f2c45f1SShri Abhyankar 
7895f2c45f1SShri Abhyankar   Level: intermediate
7905f2c45f1SShri Abhyankar 
7915f2c45f1SShri Abhyankar .seealso: DMNetworkCreate
7925f2c45f1SShri Abhyankar @*/
793d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap)
7945f2c45f1SShri Abhyankar {
795d3464fd4SAdrian Maldonado   MPI_Comm       comm;
7965f2c45f1SShri Abhyankar   PetscErrorCode ierr;
797d3464fd4SAdrian Maldonado   PetscMPIInt    size;
798d3464fd4SAdrian Maldonado   DM_Network     *oldDMnetwork = (DM_Network*)((*dm)->data);
799d3464fd4SAdrian Maldonado   DM_Network     *newDMnetwork;
8005f2c45f1SShri Abhyankar   PetscSF        pointsf;
8015f2c45f1SShri Abhyankar   DM             newDM;
80251ac5effSHong Zhang   PetscPartitioner part;
8035f2c45f1SShri Abhyankar 
8045f2c45f1SShri Abhyankar   PetscFunctionBegin;
805d3464fd4SAdrian Maldonado 
806d3464fd4SAdrian Maldonado   ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr);
807d3464fd4SAdrian Maldonado   ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr);
808d3464fd4SAdrian Maldonado   if (size == 1) PetscFunctionReturn(0);
809d3464fd4SAdrian Maldonado 
810d3464fd4SAdrian Maldonado   ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr);
8115f2c45f1SShri Abhyankar   newDMnetwork = (DM_Network*)newDM->data;
8125f2c45f1SShri Abhyankar   newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType);
81351ac5effSHong Zhang 
81451ac5effSHong Zhang   /* Enable runtime options for petscpartitioner */
81551ac5effSHong Zhang   ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr);
81651ac5effSHong Zhang   ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr);
81751ac5effSHong Zhang 
8185f2c45f1SShri Abhyankar   /* Distribute plex dm and dof section */
81980cf41d5SMatthew G. Knepley   ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr);
82051ac5effSHong Zhang 
8215f2c45f1SShri Abhyankar   /* Distribute dof section */
822d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr);
8235f2c45f1SShri Abhyankar   ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr);
824d3464fd4SAdrian Maldonado   ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr);
82551ac5effSHong Zhang 
8265f2c45f1SShri Abhyankar   /* Distribute data and associated section */
82731da1fc8SHong Zhang   ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr);
82824121865SAdrian Maldonado 
8295f2c45f1SShri Abhyankar   ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr);
8305f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr);
8315f2c45f1SShri Abhyankar   ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr);
8325f2c45f1SShri Abhyankar   newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart;
833*6fefedf4SHong Zhang   newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart;
834*6fefedf4SHong Zhang   newDMnetwork->NVertices = oldDMnetwork->NVertices;
8355f2c45f1SShri Abhyankar   newDMnetwork->NEdges = oldDMnetwork->NEdges;
83624121865SAdrian Maldonado 
8375f2c45f1SShri Abhyankar   /* Set Dof section as the default section for dm */
8385f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr);
8395f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr);
8405f2c45f1SShri Abhyankar 
84124121865SAdrian Maldonado   /* Destroy point SF */
84224121865SAdrian Maldonado   ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr);
84324121865SAdrian Maldonado 
844d3464fd4SAdrian Maldonado   ierr = DMDestroy(dm);CHKERRQ(ierr);
845d3464fd4SAdrian Maldonado   *dm  = newDM;
8465f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
8475f2c45f1SShri Abhyankar }
8485f2c45f1SShri Abhyankar 
84924121865SAdrian Maldonado /*@C
85024121865SAdrian Maldonado   PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering.
85124121865SAdrian Maldonado 
85224121865SAdrian Maldonado   Input Parameters:
85324121865SAdrian Maldonado + masterSF - the original SF structure
85424121865SAdrian Maldonado - map      - a ISLocalToGlobal mapping that contains the subset of points
85524121865SAdrian Maldonado 
85624121865SAdrian Maldonado   Output Parameters:
85724121865SAdrian Maldonado . subSF    - a subset of the masterSF for the desired subset.
85824121865SAdrian Maldonado */
85924121865SAdrian Maldonado 
86024121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) {
86124121865SAdrian Maldonado 
86224121865SAdrian Maldonado   PetscErrorCode        ierr;
86324121865SAdrian Maldonado   PetscInt              nroots, nleaves, *ilocal_sub;
86424121865SAdrian Maldonado   PetscInt              i, *ilocal_map, nroots_sub, nleaves_sub = 0;
86524121865SAdrian Maldonado   PetscInt              *local_points, *remote_points;
86624121865SAdrian Maldonado   PetscSFNode           *iremote_sub;
86724121865SAdrian Maldonado   const PetscInt        *ilocal;
86824121865SAdrian Maldonado   const PetscSFNode     *iremote;
86924121865SAdrian Maldonado 
87024121865SAdrian Maldonado   PetscFunctionBegin;
87124121865SAdrian Maldonado   ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr);
87224121865SAdrian Maldonado 
87324121865SAdrian Maldonado   /* Look for leaves that pertain to the subset of points. Get the local ordering */
87424121865SAdrian Maldonado   ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr);
87524121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr);
87624121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
87724121865SAdrian Maldonado     if (ilocal_map[i] != -1) nleaves_sub += 1;
87824121865SAdrian Maldonado   }
87924121865SAdrian Maldonado   /* Re-number ilocal with subset numbering. Need information from roots */
88024121865SAdrian Maldonado   ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr);
88124121865SAdrian Maldonado   for (i = 0; i < nroots; i++) local_points[i] = i;
88224121865SAdrian Maldonado   ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr);
88324121865SAdrian Maldonado   ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
88424121865SAdrian Maldonado   ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr);
88524121865SAdrian Maldonado   /* Fill up graph using local (that is, local to the subset) numbering. */
8864b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr);
8874b70a8deSAdrian Maldonado   ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr);
88824121865SAdrian Maldonado   nleaves_sub = 0;
88924121865SAdrian Maldonado   for (i = 0; i < nleaves; i++) {
89024121865SAdrian Maldonado     if (ilocal_map[i] != -1) {
89124121865SAdrian Maldonado       ilocal_sub[nleaves_sub] = ilocal_map[i];
8924b70a8deSAdrian Maldonado       iremote_sub[nleaves_sub].rank = iremote[i].rank;
89324121865SAdrian Maldonado       iremote_sub[nleaves_sub].index = remote_points[ilocal[i]];
89424121865SAdrian Maldonado       nleaves_sub += 1;
89524121865SAdrian Maldonado     }
89624121865SAdrian Maldonado   }
89724121865SAdrian Maldonado   ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr);
89824121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr);
89924121865SAdrian Maldonado 
90024121865SAdrian Maldonado   /* Create new subSF */
90124121865SAdrian Maldonado   ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr);
90224121865SAdrian Maldonado   ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr);
9034b70a8deSAdrian Maldonado   ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr);
90424121865SAdrian Maldonado   ierr = PetscFree(ilocal_map);CHKERRQ(ierr);
9054b70a8deSAdrian Maldonado   ierr = PetscFree(iremote_sub);CHKERRQ(ierr);
90624121865SAdrian Maldonado   PetscFunctionReturn(0);
90724121865SAdrian Maldonado }
90824121865SAdrian Maldonado 
9095f2c45f1SShri Abhyankar /*@C
9105f2c45f1SShri Abhyankar   DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point
9115f2c45f1SShri Abhyankar 
9125f2c45f1SShri Abhyankar   Not Collective
9135f2c45f1SShri Abhyankar 
9145f2c45f1SShri Abhyankar   Input Parameters:
9155f2c45f1SShri Abhyankar + dm - The DMNetwork object
9165f2c45f1SShri Abhyankar - p  - the vertex point
9175f2c45f1SShri Abhyankar 
9185f2c45f1SShri Abhyankar   Output Paramters:
9195f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point
9205f2c45f1SShri Abhyankar - edges  - List of edge points
9215f2c45f1SShri Abhyankar 
9225f2c45f1SShri Abhyankar   Level: intermediate
9235f2c45f1SShri Abhyankar 
9245f2c45f1SShri Abhyankar   Fortran Notes:
9255f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
9265f2c45f1SShri Abhyankar   include petsc.h90 in your code.
9275f2c45f1SShri Abhyankar 
928d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices
9295f2c45f1SShri Abhyankar @*/
9305f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[])
9315f2c45f1SShri Abhyankar {
9325f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9335f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9345f2c45f1SShri Abhyankar 
9355f2c45f1SShri Abhyankar   PetscFunctionBegin;
9365f2c45f1SShri Abhyankar   ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr);
9375f2c45f1SShri Abhyankar   ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);
9385f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9395f2c45f1SShri Abhyankar }
9405f2c45f1SShri Abhyankar 
9415f2c45f1SShri Abhyankar /*@C
942d842c372SHong Zhang   DMNetworkGetConnectedVertices - Return the connected vertices for this edge point
9435f2c45f1SShri Abhyankar 
9445f2c45f1SShri Abhyankar   Not Collective
9455f2c45f1SShri Abhyankar 
9465f2c45f1SShri Abhyankar   Input Parameters:
9475f2c45f1SShri Abhyankar + dm - The DMNetwork object
9485f2c45f1SShri Abhyankar - p  - the edge point
9495f2c45f1SShri Abhyankar 
9505f2c45f1SShri Abhyankar   Output Paramters:
9515f2c45f1SShri Abhyankar . vertices  - vertices connected to this edge
9525f2c45f1SShri Abhyankar 
9535f2c45f1SShri Abhyankar   Level: intermediate
9545f2c45f1SShri Abhyankar 
9555f2c45f1SShri Abhyankar   Fortran Notes:
9565f2c45f1SShri Abhyankar   Since it returns an array, this routine is only available in Fortran 90, and you must
9575f2c45f1SShri Abhyankar   include petsc.h90 in your code.
9585f2c45f1SShri Abhyankar 
9595f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges
9605f2c45f1SShri Abhyankar @*/
961d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[])
9625f2c45f1SShri Abhyankar {
9635f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9645f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9655f2c45f1SShri Abhyankar 
9665f2c45f1SShri Abhyankar   PetscFunctionBegin;
9675f2c45f1SShri Abhyankar   ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr);
9685f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
9695f2c45f1SShri Abhyankar }
9705f2c45f1SShri Abhyankar 
9715f2c45f1SShri Abhyankar /*@
9725f2c45f1SShri Abhyankar   DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex
9735f2c45f1SShri Abhyankar 
9745f2c45f1SShri Abhyankar   Not Collective
9755f2c45f1SShri Abhyankar 
9765f2c45f1SShri Abhyankar   Input Parameters:
9775f2c45f1SShri Abhyankar + dm - The DMNetwork object
9785f2c45f1SShri Abhyankar . p  - the vertex point
9795f2c45f1SShri Abhyankar 
9805f2c45f1SShri Abhyankar   Output Parameter:
9815f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point
9825f2c45f1SShri Abhyankar 
9835f2c45f1SShri Abhyankar   Level: intermediate
9845f2c45f1SShri Abhyankar 
985d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange
9865f2c45f1SShri Abhyankar @*/
9875f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost)
9885f2c45f1SShri Abhyankar {
9895f2c45f1SShri Abhyankar   PetscErrorCode ierr;
9905f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*)dm->data;
9915f2c45f1SShri Abhyankar   PetscInt       offsetg;
9925f2c45f1SShri Abhyankar   PetscSection   sectiong;
9935f2c45f1SShri Abhyankar 
9945f2c45f1SShri Abhyankar   PetscFunctionBegin;
9955f2c45f1SShri Abhyankar   *isghost = PETSC_FALSE;
9965f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&sectiong);CHKERRQ(ierr);
9975f2c45f1SShri Abhyankar   ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr);
9985f2c45f1SShri Abhyankar   if (offsetg < 0) *isghost = PETSC_TRUE;
9995f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10005f2c45f1SShri Abhyankar }
10015f2c45f1SShri Abhyankar 
10025f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm)
10035f2c45f1SShri Abhyankar {
10045f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10055f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
10065f2c45f1SShri Abhyankar 
10075f2c45f1SShri Abhyankar   PetscFunctionBegin;
10085f2c45f1SShri Abhyankar   ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr);
10095f2c45f1SShri Abhyankar   ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr);
10105f2c45f1SShri Abhyankar 
10115f2c45f1SShri Abhyankar   ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr);
10125f2c45f1SShri Abhyankar   ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr);
10135f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
10145f2c45f1SShri Abhyankar }
10155f2c45f1SShri Abhyankar 
10161ad426b7SHong Zhang /*@
101717df6e9eSHong Zhang     DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices
10181ad426b7SHong Zhang                             -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)?
10191ad426b7SHong Zhang 
10201ad426b7SHong Zhang     Collective
10211ad426b7SHong Zhang 
10221ad426b7SHong Zhang     Input Parameters:
102383b2e829SHong Zhang +   dm - The DMNetwork object
102483b2e829SHong Zhang .   eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges
102583b2e829SHong Zhang -   vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices
10261ad426b7SHong Zhang 
10271ad426b7SHong Zhang     Level: intermediate
10281ad426b7SHong Zhang 
10291ad426b7SHong Zhang @*/
103083b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg)
10311ad426b7SHong Zhang {
10321ad426b7SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
10331ad426b7SHong Zhang 
10341ad426b7SHong Zhang   PetscFunctionBegin;
103583b2e829SHong Zhang   network->userEdgeJacobian   = eflg;
103683b2e829SHong Zhang   network->userVertexJacobian = vflg;
10371ad426b7SHong Zhang   PetscFunctionReturn(0);
10381ad426b7SHong Zhang }
10391ad426b7SHong Zhang 
10401ad426b7SHong Zhang /*@
104183b2e829SHong Zhang     DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network
104283b2e829SHong Zhang 
104383b2e829SHong Zhang     Not Collective
104483b2e829SHong Zhang 
104583b2e829SHong Zhang     Input Parameters:
104683b2e829SHong Zhang +   dm - The DMNetwork object
104783b2e829SHong Zhang .   p  - the edge point
10483e97b6e8SHong Zhang -   J - array (size = 3) of Jacobian submatrices for this edge point:
10493e97b6e8SHong Zhang         J[0]: this edge
1050d842c372SHong Zhang         J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices()
105183b2e829SHong Zhang 
105283b2e829SHong Zhang     Level: intermediate
105383b2e829SHong Zhang 
105483b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix
105583b2e829SHong Zhang @*/
105683b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[])
105783b2e829SHong Zhang {
105883b2e829SHong Zhang   PetscErrorCode ierr;
105983b2e829SHong Zhang   DM_Network     *network=(DM_Network*)dm->data;
106083b2e829SHong Zhang 
106183b2e829SHong Zhang   PetscFunctionBegin;
1062883e35e8SHong Zhang   if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix");
106383b2e829SHong Zhang   if (!network->Je) {
1064883e35e8SHong Zhang     ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr);
106583b2e829SHong Zhang   }
1066883e35e8SHong Zhang   network->Je[3*p]   = J[0];
1067883e35e8SHong Zhang   network->Je[3*p+1] = J[1];
1068883e35e8SHong Zhang   network->Je[3*p+2] = J[2];
106983b2e829SHong Zhang   PetscFunctionReturn(0);
107083b2e829SHong Zhang }
107183b2e829SHong Zhang 
107283b2e829SHong Zhang /*@
107376ddfea5SHong Zhang     DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network
10741ad426b7SHong Zhang 
10751ad426b7SHong Zhang     Not Collective
10761ad426b7SHong Zhang 
10771ad426b7SHong Zhang     Input Parameters:
10781ad426b7SHong Zhang +   dm - The DMNetwork object
10791ad426b7SHong Zhang .   p  - the vertex point
10803e97b6e8SHong Zhang -   J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point:
10813e97b6e8SHong Zhang         J[0]:       this vertex
10823e97b6e8SHong Zhang         J[1+2*i]:   i-th supporting edge
10833e97b6e8SHong Zhang         J[1+2*i+1]: i-th connected vertex
10841ad426b7SHong Zhang 
10851ad426b7SHong Zhang     Level: intermediate
10861ad426b7SHong Zhang 
108783b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix
10881ad426b7SHong Zhang @*/
1089883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[])
10905f2c45f1SShri Abhyankar {
10915f2c45f1SShri Abhyankar   PetscErrorCode ierr;
10925f2c45f1SShri Abhyankar   DM_Network     *network=(DM_Network*)dm->data;
1093f4431b8cSHong Zhang   PetscInt       i,*vptr,nedges,vStart,vEnd;
1094883e35e8SHong Zhang   const PetscInt *edges;
10955f2c45f1SShri Abhyankar 
10965f2c45f1SShri Abhyankar   PetscFunctionBegin;
1097883e35e8SHong Zhang   if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix");
1098883e35e8SHong Zhang 
1099883e35e8SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1100f4431b8cSHong Zhang 
110183b2e829SHong Zhang   if (!network->Jv) {
1102*6fefedf4SHong Zhang     PetscInt nVertices = network->nVertices,nedges_total;
1103883e35e8SHong Zhang     /* count nvertex_total */
11043e97b6e8SHong Zhang     nedges_total = 0;
1105883e35e8SHong Zhang     ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1106*6fefedf4SHong Zhang     ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr);
1107f4431b8cSHong Zhang 
1108883e35e8SHong Zhang     vptr[0] = 0;
1109*6fefedf4SHong Zhang     for (i=0; i<nVertices; i++) {
1110f4431b8cSHong Zhang       ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr);
1111883e35e8SHong Zhang       nedges_total += nedges;
1112f4431b8cSHong Zhang       vptr[i+1] = vptr[i] + 2*nedges + 1;
11131ad426b7SHong Zhang     }
11143e97b6e8SHong Zhang 
1115*6fefedf4SHong Zhang     ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr);
1116883e35e8SHong Zhang     network->Jvptr = vptr;
1117883e35e8SHong Zhang   }
1118883e35e8SHong Zhang 
1119883e35e8SHong Zhang   vptr = network->Jvptr;
11203e97b6e8SHong Zhang   network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */
11213e97b6e8SHong Zhang 
11223e97b6e8SHong Zhang   /* Set Jacobian for each supporting edge and connected vertex */
1123883e35e8SHong Zhang   ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr);
1124883e35e8SHong Zhang   for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i];
1125883e35e8SHong Zhang   PetscFunctionReturn(0);
1126883e35e8SHong Zhang }
1127883e35e8SHong Zhang 
1128e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
11295cf7da58SHong Zhang {
11305cf7da58SHong Zhang   PetscErrorCode ierr;
11315cf7da58SHong Zhang   PetscInt       j;
11325cf7da58SHong Zhang   PetscScalar    val=(PetscScalar)ncols;
11335cf7da58SHong Zhang 
11345cf7da58SHong Zhang   PetscFunctionBegin;
11355cf7da58SHong Zhang   if (!ghost) {
11365cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11375cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11385cf7da58SHong Zhang     }
11395cf7da58SHong Zhang   } else {
11405cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11415cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11425cf7da58SHong Zhang     }
11435cf7da58SHong Zhang   }
11445cf7da58SHong Zhang   PetscFunctionReturn(0);
11455cf7da58SHong Zhang }
11465cf7da58SHong Zhang 
1147e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
11485cf7da58SHong Zhang {
11495cf7da58SHong Zhang   PetscErrorCode ierr;
11505cf7da58SHong Zhang   PetscInt       j,ncols_u;
11515cf7da58SHong Zhang   PetscScalar    val;
11525cf7da58SHong Zhang 
11535cf7da58SHong Zhang   PetscFunctionBegin;
11545cf7da58SHong Zhang   if (!ghost) {
11555cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11565cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
11575cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
11585cf7da58SHong Zhang       ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11595cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
11605cf7da58SHong Zhang     }
11615cf7da58SHong Zhang   } else {
11625cf7da58SHong Zhang     for (j=0; j<nrows; j++) {
11635cf7da58SHong Zhang       ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
11645cf7da58SHong Zhang       val = (PetscScalar)ncols_u;
11655cf7da58SHong Zhang       ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr);
11665cf7da58SHong Zhang       ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr);
11675cf7da58SHong Zhang     }
11685cf7da58SHong Zhang   }
11695cf7da58SHong Zhang   PetscFunctionReturn(0);
11705cf7da58SHong Zhang }
11715cf7da58SHong Zhang 
1172e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz)
11735cf7da58SHong Zhang {
11745cf7da58SHong Zhang   PetscErrorCode ierr;
11755cf7da58SHong Zhang 
11765cf7da58SHong Zhang   PetscFunctionBegin;
11775cf7da58SHong Zhang   if (Ju) {
11785cf7da58SHong Zhang     ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
11795cf7da58SHong Zhang   } else {
11805cf7da58SHong Zhang     ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr);
11815cf7da58SHong Zhang   }
11825cf7da58SHong Zhang   PetscFunctionReturn(0);
11835cf7da58SHong Zhang }
11845cf7da58SHong Zhang 
1185e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1186883e35e8SHong Zhang {
1187883e35e8SHong Zhang   PetscErrorCode ierr;
1188883e35e8SHong Zhang   PetscInt       j,*cols;
1189883e35e8SHong Zhang   PetscScalar    *zeros;
1190883e35e8SHong Zhang 
1191883e35e8SHong Zhang   PetscFunctionBegin;
1192883e35e8SHong Zhang   ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr);
1193883e35e8SHong Zhang   for (j=0; j<ncols; j++) cols[j] = j+ cstart;
1194883e35e8SHong Zhang   ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr);
1195883e35e8SHong Zhang   ierr = PetscFree2(cols,zeros);CHKERRQ(ierr);
11961ad426b7SHong Zhang   PetscFunctionReturn(0);
11971ad426b7SHong Zhang }
1198a4e85ca8SHong Zhang 
1199e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
12003e97b6e8SHong Zhang {
12013e97b6e8SHong Zhang   PetscErrorCode ierr;
12023e97b6e8SHong Zhang   PetscInt       j,M,N,row,col,ncols_u;
12033e97b6e8SHong Zhang   const PetscInt *cols;
12043e97b6e8SHong Zhang   PetscScalar    zero=0.0;
12053e97b6e8SHong Zhang 
12063e97b6e8SHong Zhang   PetscFunctionBegin;
12073e97b6e8SHong Zhang   ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr);
12083e97b6e8SHong 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);
12093e97b6e8SHong Zhang 
12103e97b6e8SHong Zhang   for (row=0; row<nrows; row++) {
12113e97b6e8SHong Zhang     ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
12123e97b6e8SHong Zhang     for (j=0; j<ncols_u; j++) {
12133e97b6e8SHong Zhang       col = cols[j] + cstart;
12143e97b6e8SHong Zhang       ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr);
12153e97b6e8SHong Zhang     }
12163e97b6e8SHong Zhang     ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr);
12173e97b6e8SHong Zhang   }
12183e97b6e8SHong Zhang   PetscFunctionReturn(0);
12193e97b6e8SHong Zhang }
12201ad426b7SHong Zhang 
1221e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J)
1222a4e85ca8SHong Zhang {
1223a4e85ca8SHong Zhang   PetscErrorCode ierr;
1224f4431b8cSHong Zhang 
1225a4e85ca8SHong Zhang   PetscFunctionBegin;
1226a4e85ca8SHong Zhang   if (Ju) {
1227a4e85ca8SHong Zhang     ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1228a4e85ca8SHong Zhang   } else {
1229a4e85ca8SHong Zhang     ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1230a4e85ca8SHong Zhang   }
1231a4e85ca8SHong Zhang   PetscFunctionReturn(0);
1232a4e85ca8SHong Zhang }
1233a4e85ca8SHong Zhang 
123424121865SAdrian 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.
123524121865SAdrian Maldonado */
123624121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog)
123724121865SAdrian Maldonado {
123824121865SAdrian Maldonado   PetscErrorCode ierr;
123924121865SAdrian Maldonado   PetscInt       i, size, dof;
124024121865SAdrian Maldonado   PetscInt       *glob2loc;
124124121865SAdrian Maldonado 
124224121865SAdrian Maldonado   PetscFunctionBegin;
124324121865SAdrian Maldonado   ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr);
124424121865SAdrian Maldonado   ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr);
124524121865SAdrian Maldonado 
124624121865SAdrian Maldonado   for (i = 0; i < size; i++) {
124724121865SAdrian Maldonado     ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr);
124824121865SAdrian Maldonado     dof = (dof >= 0) ? dof : -(dof + 1);
124924121865SAdrian Maldonado     glob2loc[i] = dof;
125024121865SAdrian Maldonado   }
125124121865SAdrian Maldonado 
125224121865SAdrian Maldonado   ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr);
125324121865SAdrian Maldonado #if 0
125424121865SAdrian Maldonado   ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr);
125524121865SAdrian Maldonado #endif
125624121865SAdrian Maldonado   PetscFunctionReturn(0);
125724121865SAdrian Maldonado }
125824121865SAdrian Maldonado 
12591ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J)
12601ad426b7SHong Zhang {
12611ad426b7SHong Zhang   PetscErrorCode ierr;
126224121865SAdrian Maldonado   PetscMPIInt    rank, size;
12631ad426b7SHong Zhang   DM_Network     *network = (DM_Network*) dm->data;
1264a4e85ca8SHong Zhang   PetscInt       eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize;
1265840c2264SHong Zhang   PetscInt       cstart,ncols,j,e,v;
126624121865SAdrian Maldonado   PetscBool      ghost,ghost_vc,ghost2,isNest;
1267a4e85ca8SHong Zhang   Mat            Juser;
1268bfbc38dcSHong Zhang   PetscSection   sectionGlobal;
1269447d78afSSatish Balay   PetscInt       nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */
1270a4e85ca8SHong Zhang   const PetscInt *edges,*cone;
12715cf7da58SHong Zhang   MPI_Comm       comm;
127224121865SAdrian Maldonado   MatType        mtype;
12735cf7da58SHong Zhang   Vec            vd_nz,vo_nz;
12745cf7da58SHong Zhang   PetscInt       *dnnz,*onnz;
12755cf7da58SHong Zhang   PetscScalar    *vdnz,*vonz;
12761ad426b7SHong Zhang 
12771ad426b7SHong Zhang   PetscFunctionBegin;
127824121865SAdrian Maldonado   mtype = dm->mattype;
127924121865SAdrian Maldonado   ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr);
128024121865SAdrian Maldonado 
128124121865SAdrian Maldonado   if (isNest) {
12820731d606SHong Zhang     /* ierr = DMCreateMatrix_Network_Nest(); */
128324121865SAdrian Maldonado     PetscInt   eDof, vDof;
128424121865SAdrian Maldonado     Mat        j11, j12, j21, j22, bA[2][2];
128524121865SAdrian Maldonado     ISLocalToGlobalMapping eISMap, vISMap;
128624121865SAdrian Maldonado 
128724121865SAdrian Maldonado     ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
128824121865SAdrian Maldonado     ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr);
128924121865SAdrian Maldonado     ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr);
129024121865SAdrian Maldonado 
129124121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr);
129224121865SAdrian Maldonado     ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr);
129324121865SAdrian Maldonado 
129424121865SAdrian Maldonado     ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr);
129524121865SAdrian Maldonado     ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
129624121865SAdrian Maldonado     ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr);
129724121865SAdrian Maldonado 
129824121865SAdrian Maldonado     ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr);
129924121865SAdrian Maldonado     ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr);
130024121865SAdrian Maldonado     ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr);
130124121865SAdrian Maldonado 
130224121865SAdrian Maldonado     ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr);
130324121865SAdrian Maldonado     ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
130424121865SAdrian Maldonado     ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr);
130524121865SAdrian Maldonado 
130624121865SAdrian Maldonado     ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr);
130724121865SAdrian Maldonado     ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr);
130824121865SAdrian Maldonado     ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr);
130924121865SAdrian Maldonado 
131024121865SAdrian Maldonado     bA[0][0] = j11;
131124121865SAdrian Maldonado     bA[0][1] = j12;
131224121865SAdrian Maldonado     bA[1][0] = j21;
131324121865SAdrian Maldonado     bA[1][1] = j22;
131424121865SAdrian Maldonado 
131524121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr);
131624121865SAdrian Maldonado     ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr);
131724121865SAdrian Maldonado 
131824121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr);
131924121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr);
132024121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr);
132124121865SAdrian Maldonado     ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr);
132224121865SAdrian Maldonado 
132324121865SAdrian Maldonado     ierr = MatSetUp(j11);CHKERRQ(ierr);
132424121865SAdrian Maldonado     ierr = MatSetUp(j12);CHKERRQ(ierr);
132524121865SAdrian Maldonado     ierr = MatSetUp(j21);CHKERRQ(ierr);
132624121865SAdrian Maldonado     ierr = MatSetUp(j22);CHKERRQ(ierr);
132724121865SAdrian Maldonado 
132824121865SAdrian Maldonado     ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr);
132924121865SAdrian Maldonado     ierr = MatSetUp(*J);CHKERRQ(ierr);
133024121865SAdrian Maldonado     ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr);
133124121865SAdrian Maldonado     ierr = MatDestroy(&j11);CHKERRQ(ierr);
133224121865SAdrian Maldonado     ierr = MatDestroy(&j12);CHKERRQ(ierr);
133324121865SAdrian Maldonado     ierr = MatDestroy(&j21);CHKERRQ(ierr);
133424121865SAdrian Maldonado     ierr = MatDestroy(&j22);CHKERRQ(ierr);
133524121865SAdrian Maldonado 
133624121865SAdrian Maldonado     ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
133724121865SAdrian Maldonado     ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
133824121865SAdrian Maldonado 
133924121865SAdrian Maldonado     /* Free structures */
134024121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr);
134124121865SAdrian Maldonado     ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr);
134224121865SAdrian Maldonado 
134324121865SAdrian Maldonado     PetscFunctionReturn(0);
134424121865SAdrian Maldonado   } else if (!network->userEdgeJacobian && !network->userVertexJacobian) {
1345a4e85ca8SHong Zhang     /* user does not provide Jacobian blocks */
1346bfbc38dcSHong Zhang     ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr);
1347bfbc38dcSHong Zhang     ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
13481ad426b7SHong Zhang     PetscFunctionReturn(0);
13491ad426b7SHong Zhang   }
13501ad426b7SHong Zhang 
1351bfbc38dcSHong Zhang   ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr);
13522a945128SHong Zhang   ierr = DMGetDefaultGlobalSection(network->plex,&sectionGlobal);CHKERRQ(ierr);
1353bfbc38dcSHong Zhang   ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr);
1354bfbc38dcSHong Zhang   ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr);
13552a945128SHong Zhang 
13562a945128SHong Zhang   ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr);
13572a945128SHong Zhang   ierr = MatSetFromOptions(*J);CHKERRQ(ierr);
135889898e50SHong Zhang 
135989898e50SHong Zhang   /* (1) Set matrix preallocation */
136089898e50SHong Zhang   /*------------------------------*/
1361840c2264SHong Zhang   ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr);
1362840c2264SHong Zhang   ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr);
1363840c2264SHong Zhang   ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr);
1364840c2264SHong Zhang   ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr);
1365840c2264SHong Zhang   ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr);
1366840c2264SHong Zhang   ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr);
1367840c2264SHong Zhang 
136889898e50SHong Zhang   /* Set preallocation for edges */
136989898e50SHong Zhang   /*-----------------------------*/
1370840c2264SHong Zhang   ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr);
1371840c2264SHong Zhang 
1372bdcb62a2SHong Zhang   ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr);
1373840c2264SHong Zhang   for (e=eStart; e<eEnd; e++) {
1374840c2264SHong Zhang     /* Get row indices */
1375840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
1376840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
1377840c2264SHong Zhang     if (nrows) {
1378840c2264SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1379840c2264SHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
1380840c2264SHong Zhang 
13815cf7da58SHong Zhang       /* Set preallocation for conntected vertices */
1382d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1383840c2264SHong Zhang       for (v=0; v<2; v++) {
1384840c2264SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
1385840c2264SHong Zhang 
1386840c2264SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1387840c2264SHong Zhang         ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr);
13885cf7da58SHong Zhang         ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1389840c2264SHong Zhang       }
1390840c2264SHong Zhang 
139189898e50SHong Zhang       /* Set preallocation for edge self */
1392840c2264SHong Zhang       cstart = rstart;
1393840c2264SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
13945cf7da58SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1395840c2264SHong Zhang     }
1396840c2264SHong Zhang   }
1397840c2264SHong Zhang 
139889898e50SHong Zhang   /* Set preallocation for vertices */
139989898e50SHong Zhang   /*--------------------------------*/
1400840c2264SHong Zhang   ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr);
1401840c2264SHong Zhang   if (vEnd - vStart) {
1402840c2264SHong Zhang     if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv");
1403840c2264SHong Zhang     vptr = network->Jvptr;
1404840c2264SHong Zhang     if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr");
1405840c2264SHong Zhang   }
1406840c2264SHong Zhang 
1407840c2264SHong Zhang   for (v=vStart; v<vEnd; v++) {
1408840c2264SHong Zhang     /* Get row indices */
1409840c2264SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1410840c2264SHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
1411840c2264SHong Zhang     if (!nrows) continue;
1412840c2264SHong Zhang 
1413bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1414bdcb62a2SHong Zhang     if (ghost) {
1415bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1416bdcb62a2SHong Zhang     } else {
1417bdcb62a2SHong Zhang       rows_v = rows;
1418bdcb62a2SHong Zhang     }
1419bdcb62a2SHong Zhang 
1420bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1421840c2264SHong Zhang 
1422840c2264SHong Zhang     /* Get supporting edges and connected vertices */
1423840c2264SHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1424840c2264SHong Zhang 
1425840c2264SHong Zhang     for (e=0; e<nedges; e++) {
1426840c2264SHong Zhang       /* Supporting edges */
1427840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1428840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1429840c2264SHong Zhang 
1430840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1431bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr);
1432840c2264SHong Zhang 
1433840c2264SHong Zhang       /* Connected vertices */
1434d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
1435840c2264SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
1436840c2264SHong Zhang       ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr);
1437840c2264SHong Zhang 
1438840c2264SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1439840c2264SHong Zhang 
1440840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1441e102a522SHong Zhang       if (ghost_vc||ghost) {
1442e102a522SHong Zhang         ghost2 = PETSC_TRUE;
1443e102a522SHong Zhang       } else {
1444e102a522SHong Zhang         ghost2 = PETSC_FALSE;
1445e102a522SHong Zhang       }
1446e102a522SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr);
1447840c2264SHong Zhang     }
1448840c2264SHong Zhang 
144989898e50SHong Zhang     /* Set preallocation for vertex self */
1450840c2264SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1451840c2264SHong Zhang     if (!ghost) {
1452840c2264SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1453840c2264SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1454bdcb62a2SHong Zhang       ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr);
1455840c2264SHong Zhang     }
1456bdcb62a2SHong Zhang     if (ghost) {
1457bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1458bdcb62a2SHong Zhang     }
1459840c2264SHong Zhang   }
1460840c2264SHong Zhang 
1461840c2264SHong Zhang   ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr);
1462840c2264SHong Zhang   ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr);
14635cf7da58SHong Zhang 
14645cf7da58SHong Zhang   ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr);
14655cf7da58SHong Zhang 
14665cf7da58SHong Zhang   ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr);
1467840c2264SHong Zhang   ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr);
1468840c2264SHong Zhang 
1469840c2264SHong Zhang   ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr);
1470840c2264SHong Zhang   ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr);
1471840c2264SHong Zhang   for (j=0; j<localSize; j++) {
1472e102a522SHong Zhang     dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]);
1473e102a522SHong Zhang     onnz[j] = (PetscInt)PetscRealPart(vonz[j]);
1474840c2264SHong Zhang   }
1475840c2264SHong Zhang   ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr);
1476840c2264SHong Zhang   ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr);
1477840c2264SHong Zhang   ierr = VecDestroy(&vd_nz);CHKERRQ(ierr);
1478840c2264SHong Zhang   ierr = VecDestroy(&vo_nz);CHKERRQ(ierr);
1479840c2264SHong Zhang 
14805cf7da58SHong Zhang   ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr);
14815cf7da58SHong Zhang   ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr);
14825cf7da58SHong Zhang   ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr);
14835cf7da58SHong Zhang 
14845cf7da58SHong Zhang   ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr);
14855cf7da58SHong Zhang 
148689898e50SHong Zhang   /* (2) Set matrix entries for edges */
148789898e50SHong Zhang   /*----------------------------------*/
14881ad426b7SHong Zhang   for (e=eStart; e<eEnd; e++) {
1489bfbc38dcSHong Zhang     /* Get row indices */
14901ad426b7SHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr);
149117df6e9eSHong Zhang     ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr);
14924b976069SHong Zhang     if (nrows) {
14934b976069SHong Zhang       if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je");
1494bdcb62a2SHong Zhang 
149517df6e9eSHong Zhang       for (j=0; j<nrows; j++) rows[j] = j + rstart;
14961ad426b7SHong Zhang 
1497bfbc38dcSHong Zhang       /* Set matrix entries for conntected vertices */
1498d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr);
1499bfbc38dcSHong Zhang       for (v=0; v<2; v++) {
1500bfbc38dcSHong Zhang         ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr);
1501883e35e8SHong Zhang         ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr);
15023e97b6e8SHong Zhang 
1503a4e85ca8SHong Zhang         Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */
1504a4e85ca8SHong Zhang         ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr);
1505bfbc38dcSHong Zhang       }
150617df6e9eSHong Zhang 
1507bfbc38dcSHong Zhang       /* Set matrix entries for edge self */
15083e97b6e8SHong Zhang       cstart = rstart;
1509a4e85ca8SHong Zhang       Juser = network->Je[3*e]; /* Jacobian(e,e) */
1510a4e85ca8SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr);
15111ad426b7SHong Zhang     }
15124b976069SHong Zhang   }
15131ad426b7SHong Zhang 
1514bfbc38dcSHong Zhang   /* Set matrix entries for vertices */
151583b2e829SHong Zhang   /*---------------------------------*/
15161ad426b7SHong Zhang   for (v=vStart; v<vEnd; v++) {
1517bfbc38dcSHong Zhang     /* Get row indices */
1518596e729fSHong Zhang     ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr);
1519596e729fSHong Zhang     ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr);
15204b976069SHong Zhang     if (!nrows) continue;
1521596e729fSHong Zhang 
1522bdcb62a2SHong Zhang     ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr);
1523bdcb62a2SHong Zhang     if (ghost) {
1524bdcb62a2SHong Zhang       ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr);
1525bdcb62a2SHong Zhang     } else {
1526bdcb62a2SHong Zhang       rows_v = rows;
1527bdcb62a2SHong Zhang     }
1528bdcb62a2SHong Zhang     for (j=0; j<nrows; j++) rows_v[j] = j + rstart;
1529596e729fSHong Zhang 
1530bfbc38dcSHong Zhang     /* Get supporting edges and connected vertices */
1531596e729fSHong Zhang     ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr);
1532596e729fSHong Zhang 
1533596e729fSHong Zhang     for (e=0; e<nedges; e++) {
1534bfbc38dcSHong Zhang       /* Supporting edges */
1535596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr);
1536596e729fSHong Zhang       ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr);
1537596e729fSHong Zhang 
1538a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */
1539bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1540596e729fSHong Zhang 
1541bfbc38dcSHong Zhang       /* Connected vertices */
1542d842c372SHong Zhang       ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr);
15432a945128SHong Zhang       vc = (v == cone[0]) ? cone[1]:cone[0];
15442a945128SHong Zhang 
154544aca652SHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr);
154644aca652SHong Zhang       ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr);
1547a4e85ca8SHong Zhang 
1548a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */
1549bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr);
1550596e729fSHong Zhang     }
1551596e729fSHong Zhang 
1552bfbc38dcSHong Zhang     /* Set matrix entries for vertex self */
15531ad426b7SHong Zhang     if (!ghost) {
1554596e729fSHong Zhang       ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr);
1555a4e85ca8SHong Zhang       Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */
1556bdcb62a2SHong Zhang       ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr);
1557bdcb62a2SHong Zhang     }
1558bdcb62a2SHong Zhang     if (ghost) {
1559bdcb62a2SHong Zhang       ierr = PetscFree(rows_v);CHKERRQ(ierr);
1560bdcb62a2SHong Zhang     }
15611ad426b7SHong Zhang   }
1562a4e85ca8SHong Zhang   ierr = PetscFree(rows);CHKERRQ(ierr);
1563bdcb62a2SHong Zhang 
15641ad426b7SHong Zhang   ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
15651ad426b7SHong Zhang   ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr);
1566dd6f46cdSHong Zhang 
15675f2c45f1SShri Abhyankar   ierr = MatSetDM(*J,dm);CHKERRQ(ierr);
15685f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
15695f2c45f1SShri Abhyankar }
15705f2c45f1SShri Abhyankar 
15715f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm)
15725f2c45f1SShri Abhyankar {
15735f2c45f1SShri Abhyankar   PetscErrorCode ierr;
15745f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
15755f2c45f1SShri Abhyankar 
15765f2c45f1SShri Abhyankar   PetscFunctionBegin;
15778415c774SShri Abhyankar   if (--network->refct > 0) PetscFunctionReturn(0);
157883b2e829SHong Zhang   if (network->Je) {
157983b2e829SHong Zhang     ierr = PetscFree(network->Je);CHKERRQ(ierr);
158083b2e829SHong Zhang   }
158183b2e829SHong Zhang   if (network->Jv) {
1582883e35e8SHong Zhang     ierr = PetscFree(network->Jvptr);CHKERRQ(ierr);
158383b2e829SHong Zhang     ierr = PetscFree(network->Jv);CHKERRQ(ierr);
15841ad426b7SHong Zhang   }
158513c2a604SAdrian Maldonado 
158613c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr);
158713c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr);
158813c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr);
158913c2a604SAdrian Maldonado   if (network->vertex.sf) {
159013c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr);
159113c2a604SAdrian Maldonado   }
159213c2a604SAdrian Maldonado   /* edge */
159313c2a604SAdrian Maldonado   ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr);
159413c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr);
159513c2a604SAdrian Maldonado   ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr);
159613c2a604SAdrian Maldonado   if (network->edge.sf) {
159713c2a604SAdrian Maldonado     ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr);
159813c2a604SAdrian Maldonado   }
15995f2c45f1SShri Abhyankar   ierr = DMDestroy(&network->plex);CHKERRQ(ierr);
16005f2c45f1SShri Abhyankar   network->edges = NULL;
16015f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr);
16025f2c45f1SShri Abhyankar   ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr);
160383b2e829SHong Zhang 
16045f2c45f1SShri Abhyankar   ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr);
16055f2c45f1SShri Abhyankar   ierr = PetscFree(network->cvalue);CHKERRQ(ierr);
16065f2c45f1SShri Abhyankar   ierr = PetscFree(network->header);CHKERRQ(ierr);
16075f2c45f1SShri Abhyankar   ierr = PetscFree(network);CHKERRQ(ierr);
16085f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16095f2c45f1SShri Abhyankar }
16105f2c45f1SShri Abhyankar 
16115f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer)
16125f2c45f1SShri Abhyankar {
16135f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16145f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16155f2c45f1SShri Abhyankar 
16165f2c45f1SShri Abhyankar   PetscFunctionBegin;
16175f2c45f1SShri Abhyankar   ierr = DMView(network->plex,viewer);CHKERRQ(ierr);
16185f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16195f2c45f1SShri Abhyankar }
16205f2c45f1SShri Abhyankar 
16215f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l)
16225f2c45f1SShri Abhyankar {
16235f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16245f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16255f2c45f1SShri Abhyankar 
16265f2c45f1SShri Abhyankar   PetscFunctionBegin;
16275f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr);
16285f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16295f2c45f1SShri Abhyankar }
16305f2c45f1SShri Abhyankar 
16315f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l)
16325f2c45f1SShri Abhyankar {
16335f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16345f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16355f2c45f1SShri Abhyankar 
16365f2c45f1SShri Abhyankar   PetscFunctionBegin;
16375f2c45f1SShri Abhyankar   ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr);
16385f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16395f2c45f1SShri Abhyankar }
16405f2c45f1SShri Abhyankar 
16415f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g)
16425f2c45f1SShri Abhyankar {
16435f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16445f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16455f2c45f1SShri Abhyankar 
16465f2c45f1SShri Abhyankar   PetscFunctionBegin;
16475f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr);
16485f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16495f2c45f1SShri Abhyankar }
16505f2c45f1SShri Abhyankar 
16515f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g)
16525f2c45f1SShri Abhyankar {
16535f2c45f1SShri Abhyankar   PetscErrorCode ierr;
16545f2c45f1SShri Abhyankar   DM_Network     *network = (DM_Network*) dm->data;
16555f2c45f1SShri Abhyankar 
16565f2c45f1SShri Abhyankar   PetscFunctionBegin;
16575f2c45f1SShri Abhyankar   ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr);
16585f2c45f1SShri Abhyankar   PetscFunctionReturn(0);
16595f2c45f1SShri Abhyankar }
1660