15f2c45f1SShri Abhyankar #include <petsc-private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 25f2c45f1SShri Abhyankar #include <petscdmplex.h> 35f2c45f1SShri Abhyankar #include <petscsf.h> 45f2c45f1SShri Abhyankar 55f2c45f1SShri Abhyankar #undef __FUNCT__ 65f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkSetSizes" 75f2c45f1SShri Abhyankar /*@ 85f2c45f1SShri Abhyankar DMNetworkSetSizes - Sets the local and global vertices and edges. 95f2c45f1SShri Abhyankar 105f2c45f1SShri Abhyankar Collective on DM 115f2c45f1SShri Abhyankar 125f2c45f1SShri Abhyankar Input Parameters: 135f2c45f1SShri Abhyankar + dm - the dm object 145f2c45f1SShri Abhyankar . nV - number of local vertices 155f2c45f1SShri Abhyankar . nE - number of local edges 165f2c45f1SShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) 175f2c45f1SShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE) 185f2c45f1SShri Abhyankar 195f2c45f1SShri Abhyankar Notes 205f2c45f1SShri Abhyankar If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 215f2c45f1SShri Abhyankar 225f2c45f1SShri Abhyankar You cannot change the sizes once they have been set 235f2c45f1SShri Abhyankar 245f2c45f1SShri Abhyankar Level: intermediate 255f2c45f1SShri Abhyankar 265f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 275f2c45f1SShri Abhyankar @*/ 285f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE) 295f2c45f1SShri Abhyankar { 305f2c45f1SShri Abhyankar PetscErrorCode ierr; 315f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 325f2c45f1SShri Abhyankar PetscInt a[2],b[2]; 335f2c45f1SShri Abhyankar 345f2c45f1SShri Abhyankar PetscFunctionBegin; 355f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 365f2c45f1SShri Abhyankar if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4); 375f2c45f1SShri Abhyankar if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5); 385f2c45f1SShri Abhyankar if (NV > 0 && nV > NV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local vertex size %D cannot be larger than global vertex size %D",nV,NV); 395f2c45f1SShri Abhyankar if (NE > 0 && nE > NE) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local edge size %D cannot be larger than global edge size %D",nE,NE); 405f2c45f1SShri Abhyankar if ((network->nNodes >= 0 || network->NNodes >= 0) && (network->nNodes != nV || network->NNodes != NV)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset vertex sizes to %D local %D global after previously setting them to %D local %D global",nV,NV,network->nNodes,network->NNodes); 415f2c45f1SShri Abhyankar if ((network->nEdges >= 0 || network->NEdges >= 0) && (network->nEdges != nE || network->NEdges != NE)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset edge sizes to %D local %D global after previously setting them to %D local %D global",nE,NE,network->nEdges,network->NEdges); 425f2c45f1SShri Abhyankar if (NE < 0 || NV < 0) { 435f2c45f1SShri Abhyankar a[0] = nV; a[1] = nE; 445f2c45f1SShri Abhyankar ierr = MPI_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 455f2c45f1SShri Abhyankar NV = b[0]; NE = b[1]; 465f2c45f1SShri Abhyankar } 475f2c45f1SShri Abhyankar network->nNodes = nV; 485f2c45f1SShri Abhyankar network->NNodes = NV; 495f2c45f1SShri Abhyankar network->nEdges = nE; 505f2c45f1SShri Abhyankar network->NEdges = NE; 515f2c45f1SShri Abhyankar PetscFunctionReturn(0); 525f2c45f1SShri Abhyankar } 535f2c45f1SShri Abhyankar 545f2c45f1SShri Abhyankar #undef __FUNCT__ 555f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkSetEdgeList" 565f2c45f1SShri Abhyankar /*@ 575f2c45f1SShri Abhyankar DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 585f2c45f1SShri Abhyankar 595f2c45f1SShri Abhyankar Logically collective on DM 605f2c45f1SShri Abhyankar 615f2c45f1SShri Abhyankar Input Parameters: 625f2c45f1SShri Abhyankar . edges - list of edges 635f2c45f1SShri Abhyankar 645f2c45f1SShri Abhyankar Notes: 655f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 665f2c45f1SShri Abhyankar not be destroyed before the call to DMNetworkLayoutSetUp 675f2c45f1SShri Abhyankar 685f2c45f1SShri Abhyankar Level: intermediate 695f2c45f1SShri Abhyankar 705f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes 715f2c45f1SShri Abhyankar @*/ 725f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[]) 735f2c45f1SShri Abhyankar { 745f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 755f2c45f1SShri Abhyankar 765f2c45f1SShri Abhyankar PetscFunctionBegin; 775f2c45f1SShri Abhyankar network->edges = edgelist; 785f2c45f1SShri Abhyankar PetscFunctionReturn(0); 795f2c45f1SShri Abhyankar } 805f2c45f1SShri Abhyankar 815f2c45f1SShri Abhyankar #undef __FUNCT__ 825f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkLayoutSetUp" 835f2c45f1SShri Abhyankar /*@ 845f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 855f2c45f1SShri Abhyankar 865f2c45f1SShri Abhyankar Collective on DM 875f2c45f1SShri Abhyankar 885f2c45f1SShri Abhyankar Input Parameters 895f2c45f1SShri Abhyankar . DM - the dmnetwork object 905f2c45f1SShri Abhyankar 915f2c45f1SShri Abhyankar Notes: 925f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 935f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 945f2c45f1SShri Abhyankar 955f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 965f2c45f1SShri Abhyankar 975f2c45f1SShri Abhyankar Level: intermediate 985f2c45f1SShri Abhyankar 995f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 1005f2c45f1SShri Abhyankar @*/ 1015f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 1025f2c45f1SShri Abhyankar { 1035f2c45f1SShri Abhyankar PetscErrorCode ierr; 1045f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1055f2c45f1SShri Abhyankar PetscInt dim = 1; /* One dimensional network */ 1065f2c45f1SShri Abhyankar PetscInt numCorners=2; 1075f2c45f1SShri Abhyankar PetscInt spacedim=2; 1085f2c45f1SShri Abhyankar double *vertexcoords=NULL; 1095f2c45f1SShri Abhyankar PetscInt i; 1105f2c45f1SShri Abhyankar PetscInt ndata; 1115f2c45f1SShri Abhyankar 1125f2c45f1SShri Abhyankar PetscFunctionBegin; 1135f2c45f1SShri Abhyankar if (network->nNodes) { 114*75b160a0SShri Abhyankar ierr = PetscMalloc1(numCorners*network->nNodes,&vertexcoords);CHKERRQ(ierr); 1155f2c45f1SShri Abhyankar } 1165f2c45f1SShri Abhyankar ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 1175f2c45f1SShri Abhyankar if (network->nNodes) { 1185f2c45f1SShri Abhyankar ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 1195f2c45f1SShri Abhyankar } 1205f2c45f1SShri Abhyankar ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 1215f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 1225f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 1235f2c45f1SShri Abhyankar 1245f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 1255f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 1265f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 1275f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 1285f2c45f1SShri Abhyankar 1295f2c45f1SShri Abhyankar network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 130*75b160a0SShri Abhyankar ierr = PetscMalloc1((network->pEnd-network->pStart),&network->header);CHKERRQ(ierr); 1315f2c45f1SShri Abhyankar for (i = network->pStart; i < network->pEnd; i++) { 1325f2c45f1SShri Abhyankar network->header[i].ndata = 0; 1335f2c45f1SShri Abhyankar ndata = network->header[i].ndata; 1345f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 1355f2c45f1SShri Abhyankar network->header[i].offset[ndata] = 0; 1365f2c45f1SShri Abhyankar } 137*75b160a0SShri Abhyankar ierr = PetscMalloc1((network->pEnd-network->pStart),&network->cvalue);CHKERRQ(ierr); 1385f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1395f2c45f1SShri Abhyankar } 1405f2c45f1SShri Abhyankar 1415f2c45f1SShri Abhyankar #undef __FUNCT__ 1425f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkRegisterComponent" 1435f2c45f1SShri Abhyankar /*@ 1445f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 1455f2c45f1SShri Abhyankar 1465f2c45f1SShri Abhyankar Logically collective on DM 1475f2c45f1SShri Abhyankar 1485f2c45f1SShri Abhyankar Input Parameters 1495f2c45f1SShri Abhyankar + dm - the network object 1505f2c45f1SShri Abhyankar . name - the component name 1515f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 1525f2c45f1SShri Abhyankar 1535f2c45f1SShri Abhyankar Output Parameters 1545f2c45f1SShri Abhyankar . key - an integer key that defines the component 1555f2c45f1SShri Abhyankar 1565f2c45f1SShri Abhyankar Notes 1575f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 1585f2c45f1SShri Abhyankar 1595f2c45f1SShri Abhyankar Level: intermediate 1605f2c45f1SShri Abhyankar 1615f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 1625f2c45f1SShri Abhyankar @*/ 1635f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 1645f2c45f1SShri Abhyankar { 1655f2c45f1SShri Abhyankar PetscErrorCode ierr; 1665f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1675f2c45f1SShri Abhyankar DMNetworkComponent *component=&network->component[network->ncomponent]; 1685f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 1695f2c45f1SShri Abhyankar PetscInt i; 1705f2c45f1SShri Abhyankar 1715f2c45f1SShri Abhyankar PetscFunctionBegin; 1725f2c45f1SShri Abhyankar 1735f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 1745f2c45f1SShri Abhyankar ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 1755f2c45f1SShri Abhyankar if (flg) { 1765f2c45f1SShri Abhyankar *key = i; 1775f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1785f2c45f1SShri Abhyankar } 1795f2c45f1SShri Abhyankar } 1805f2c45f1SShri Abhyankar 1815f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 1825f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 1835f2c45f1SShri Abhyankar *key = network->ncomponent; 1845f2c45f1SShri Abhyankar network->ncomponent++; 1855f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1865f2c45f1SShri Abhyankar } 1875f2c45f1SShri Abhyankar 1885f2c45f1SShri Abhyankar #undef __FUNCT__ 1895f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVertexRange" 1905f2c45f1SShri Abhyankar /*@ 1915f2c45f1SShri Abhyankar DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 1925f2c45f1SShri Abhyankar 1935f2c45f1SShri Abhyankar Not Collective 1945f2c45f1SShri Abhyankar 1955f2c45f1SShri Abhyankar Input Parameters: 1965f2c45f1SShri Abhyankar + dm - The DMNetwork object 1975f2c45f1SShri Abhyankar 1985f2c45f1SShri Abhyankar Output Paramters: 1995f2c45f1SShri Abhyankar + vStart - The first vertex point 2005f2c45f1SShri Abhyankar - vEnd - One beyond the last vertex point 2015f2c45f1SShri Abhyankar 2025f2c45f1SShri Abhyankar Level: intermediate 2035f2c45f1SShri Abhyankar 2045f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange 2055f2c45f1SShri Abhyankar @*/ 2065f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 2075f2c45f1SShri Abhyankar { 2085f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2095f2c45f1SShri Abhyankar 2105f2c45f1SShri Abhyankar PetscFunctionBegin; 2115f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 2125f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 2135f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2145f2c45f1SShri Abhyankar } 2155f2c45f1SShri Abhyankar 2165f2c45f1SShri Abhyankar #undef __FUNCT__ 2175f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetEdgeRange" 2185f2c45f1SShri Abhyankar /*@ 2195f2c45f1SShri Abhyankar DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 2205f2c45f1SShri Abhyankar 2215f2c45f1SShri Abhyankar Not Collective 2225f2c45f1SShri Abhyankar 2235f2c45f1SShri Abhyankar Input Parameters: 2245f2c45f1SShri Abhyankar + dm - The DMNetwork object 2255f2c45f1SShri Abhyankar 2265f2c45f1SShri Abhyankar Output Paramters: 2275f2c45f1SShri Abhyankar + eStart - The first edge point 2285f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 2295f2c45f1SShri Abhyankar 2305f2c45f1SShri Abhyankar Level: intermediate 2315f2c45f1SShri Abhyankar 2325f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange 2335f2c45f1SShri Abhyankar @*/ 2345f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 2355f2c45f1SShri Abhyankar { 2365f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2375f2c45f1SShri Abhyankar 2385f2c45f1SShri Abhyankar PetscFunctionBegin; 2395f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 2405f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 2415f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2425f2c45f1SShri Abhyankar } 2435f2c45f1SShri Abhyankar 2445f2c45f1SShri Abhyankar #undef __FUNCT__ 2455f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddComponent" 2465f2c45f1SShri Abhyankar /*@ 2475f2c45f1SShri Abhyankar DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 2485f2c45f1SShri Abhyankar 2495f2c45f1SShri Abhyankar Not Collective 2505f2c45f1SShri Abhyankar 2515f2c45f1SShri Abhyankar Input Parameters: 2525f2c45f1SShri Abhyankar + dm - The DMNetwork object 2535f2c45f1SShri Abhyankar . p - vertex/edge point 2545f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component 2555f2c45f1SShri Abhyankar - compvalue - pointer to the data structure for the component 2565f2c45f1SShri Abhyankar 2575f2c45f1SShri Abhyankar Level: intermediate 2585f2c45f1SShri Abhyankar 2595f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 2605f2c45f1SShri Abhyankar @*/ 2615f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 2625f2c45f1SShri Abhyankar { 2635f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2645f2c45f1SShri Abhyankar DMNetworkComponent component=network->component[componentkey]; 2655f2c45f1SShri Abhyankar DMNetworkComponentHeader header=&network->header[p]; 2665f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue=&network->cvalue[p]; 2675f2c45f1SShri Abhyankar PetscErrorCode ierr; 2685f2c45f1SShri Abhyankar 2695f2c45f1SShri Abhyankar PetscFunctionBegin; 2705f2c45f1SShri Abhyankar header->size[header->ndata] = component.size; 2715f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,p,component.size);CHKERRQ(ierr); 2725f2c45f1SShri Abhyankar header->key[header->ndata] = componentkey; 2735f2c45f1SShri Abhyankar if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 2745f2c45f1SShri Abhyankar 2755f2c45f1SShri Abhyankar cvalue->data[header->ndata] = (void*)compvalue; 2765f2c45f1SShri Abhyankar header->ndata++; 2775f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2785f2c45f1SShri Abhyankar } 2795f2c45f1SShri Abhyankar 2805f2c45f1SShri Abhyankar #undef __FUNCT__ 2815f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetNumComponents" 2825f2c45f1SShri Abhyankar /*@ 2835f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 2845f2c45f1SShri Abhyankar 2855f2c45f1SShri Abhyankar Not Collective 2865f2c45f1SShri Abhyankar 2875f2c45f1SShri Abhyankar Input Parameters: 2885f2c45f1SShri Abhyankar + dm - The DMNetwork object 2895f2c45f1SShri Abhyankar . p - vertex/edge point 2905f2c45f1SShri Abhyankar 2915f2c45f1SShri Abhyankar Output Parameters: 2925f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 2935f2c45f1SShri Abhyankar 2945f2c45f1SShri Abhyankar Level: intermediate 2955f2c45f1SShri Abhyankar 2965f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 2975f2c45f1SShri Abhyankar @*/ 2985f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 2995f2c45f1SShri Abhyankar { 3005f2c45f1SShri Abhyankar PetscErrorCode ierr; 3015f2c45f1SShri Abhyankar PetscInt offset; 3025f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3035f2c45f1SShri Abhyankar 3045f2c45f1SShri Abhyankar PetscFunctionBegin; 3055f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 3065f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 3075f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3085f2c45f1SShri Abhyankar } 3095f2c45f1SShri Abhyankar 3105f2c45f1SShri Abhyankar #undef __FUNCT__ 3115f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetComponentTypeOffset" 3125f2c45f1SShri Abhyankar /*@ 3135f2c45f1SShri Abhyankar DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the 3145f2c45f1SShri Abhyankar component value from the component data array 3155f2c45f1SShri Abhyankar 3165f2c45f1SShri Abhyankar Not Collective 3175f2c45f1SShri Abhyankar 3185f2c45f1SShri Abhyankar Input Parameters: 3195f2c45f1SShri Abhyankar + dm - The DMNetwork object 3205f2c45f1SShri Abhyankar . p - vertex/edge point 3215f2c45f1SShri Abhyankar - compnum - component number 3225f2c45f1SShri Abhyankar 3235f2c45f1SShri Abhyankar Output Parameters: 3245f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component 3255f2c45f1SShri Abhyankar - offset - offset into the component data array associated with the vertex/edge point 3265f2c45f1SShri Abhyankar 3275f2c45f1SShri Abhyankar Notes: 3285f2c45f1SShri Abhyankar Typical usage: 3295f2c45f1SShri Abhyankar 3305f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray(dm, &arr); 3315f2c45f1SShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 3325f2c45f1SShri Abhyankar Loop over vertices or edges 3335f2c45f1SShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 3345f2c45f1SShri Abhyankar Loop over numcomps 3355f2c45f1SShri Abhyankar DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset); 3365f2c45f1SShri Abhyankar compdata = (UserCompDataType)(arr+offset); 3375f2c45f1SShri Abhyankar 3385f2c45f1SShri Abhyankar Level: intermediate 3395f2c45f1SShri Abhyankar 3405f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 3415f2c45f1SShri Abhyankar @*/ 3425f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 3435f2c45f1SShri Abhyankar { 3445f2c45f1SShri Abhyankar PetscErrorCode ierr; 3455f2c45f1SShri Abhyankar PetscInt offsetp; 3465f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 3475f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3485f2c45f1SShri Abhyankar 3495f2c45f1SShri Abhyankar PetscFunctionBegin; 3505f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 3515f2c45f1SShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 3525f2c45f1SShri Abhyankar *compkey = header->key[compnum]; 3535f2c45f1SShri Abhyankar *offset = offsetp+network->dataheadersize+header->offset[compnum]; 3545f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3555f2c45f1SShri Abhyankar } 3565f2c45f1SShri Abhyankar 3575f2c45f1SShri Abhyankar #undef __FUNCT__ 3585f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVariableOffset" 3595f2c45f1SShri Abhyankar /*@ 3605f2c45f1SShri Abhyankar DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 3615f2c45f1SShri Abhyankar 3625f2c45f1SShri Abhyankar Not Collective 3635f2c45f1SShri Abhyankar 3645f2c45f1SShri Abhyankar Input Parameters: 3655f2c45f1SShri Abhyankar + dm - The DMNetwork object 3665f2c45f1SShri Abhyankar - p - the edge/vertex point 3675f2c45f1SShri Abhyankar 3685f2c45f1SShri Abhyankar Output Parameters: 3695f2c45f1SShri Abhyankar . offset - the offset 3705f2c45f1SShri Abhyankar 3715f2c45f1SShri Abhyankar Level: intermediate 3725f2c45f1SShri Abhyankar 3735f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 3745f2c45f1SShri Abhyankar @*/ 3755f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 3765f2c45f1SShri Abhyankar { 3775f2c45f1SShri Abhyankar PetscErrorCode ierr; 3785f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3795f2c45f1SShri Abhyankar 3805f2c45f1SShri Abhyankar PetscFunctionBegin; 3815f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DofSection,p,offset);CHKERRQ(ierr); 3825f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3835f2c45f1SShri Abhyankar } 3845f2c45f1SShri Abhyankar 3855f2c45f1SShri Abhyankar #undef __FUNCT__ 3865f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetVariableGlobalOffset" 3875f2c45f1SShri Abhyankar /*@ 3885f2c45f1SShri Abhyankar DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 3895f2c45f1SShri Abhyankar 3905f2c45f1SShri Abhyankar Not Collective 3915f2c45f1SShri Abhyankar 3925f2c45f1SShri Abhyankar Input Parameters: 3935f2c45f1SShri Abhyankar + dm - The DMNetwork object 3945f2c45f1SShri Abhyankar - p - the edge/vertex point 3955f2c45f1SShri Abhyankar 3965f2c45f1SShri Abhyankar Output Parameters: 3975f2c45f1SShri Abhyankar . offsetg - the offset 3985f2c45f1SShri Abhyankar 3995f2c45f1SShri Abhyankar Level: intermediate 4005f2c45f1SShri Abhyankar 4015f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 4025f2c45f1SShri Abhyankar @*/ 4035f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 4045f2c45f1SShri Abhyankar { 4055f2c45f1SShri Abhyankar PetscErrorCode ierr; 4065f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4075f2c45f1SShri Abhyankar 4085f2c45f1SShri Abhyankar PetscFunctionBegin; 4095f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);CHKERRQ(ierr); 4105f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4115f2c45f1SShri Abhyankar } 4125f2c45f1SShri Abhyankar 4135f2c45f1SShri Abhyankar #undef __FUNCT__ 4145f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddNumVariables" 4155f2c45f1SShri Abhyankar /*@ 4165f2c45f1SShri Abhyankar DMNetworkAddNumVariables - Add number of variables associated with a given point. 4175f2c45f1SShri Abhyankar 4185f2c45f1SShri Abhyankar Not Collective 4195f2c45f1SShri Abhyankar 4205f2c45f1SShri Abhyankar Input Parameters: 4215f2c45f1SShri Abhyankar + dm - The DMNetworkObject 4225f2c45f1SShri Abhyankar . p - the vertex/edge point 4235f2c45f1SShri Abhyankar - nvar - number of additional variables 4245f2c45f1SShri Abhyankar 4255f2c45f1SShri Abhyankar Level: intermediate 4265f2c45f1SShri Abhyankar 4275f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables 4285f2c45f1SShri Abhyankar @*/ 4295f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 4305f2c45f1SShri Abhyankar { 4315f2c45f1SShri Abhyankar PetscErrorCode ierr; 4325f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4335f2c45f1SShri Abhyankar 4345f2c45f1SShri Abhyankar PetscFunctionBegin; 4355f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 4365f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4375f2c45f1SShri Abhyankar } 4385f2c45f1SShri Abhyankar 4395f2c45f1SShri Abhyankar #undef __FUNCT__ 4405f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkAddNumVariables" 4415f2c45f1SShri Abhyankar /*@ 4425f2c45f1SShri Abhyankar DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 4435f2c45f1SShri Abhyankar 4445f2c45f1SShri Abhyankar Not Collective 4455f2c45f1SShri Abhyankar 4465f2c45f1SShri Abhyankar Input Parameters: 4475f2c45f1SShri Abhyankar + dm - The DMNetworkObject 4485f2c45f1SShri Abhyankar . p - the vertex/edge point 4495f2c45f1SShri Abhyankar - nvar - number of variables 4505f2c45f1SShri Abhyankar 4515f2c45f1SShri Abhyankar Level: intermediate 4525f2c45f1SShri Abhyankar 4535f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables 4545f2c45f1SShri Abhyankar @*/ 4555f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 4565f2c45f1SShri Abhyankar { 4575f2c45f1SShri Abhyankar PetscErrorCode ierr; 4585f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4595f2c45f1SShri Abhyankar 4605f2c45f1SShri Abhyankar PetscFunctionBegin; 4615f2c45f1SShri Abhyankar ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 4625f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4635f2c45f1SShri Abhyankar } 4645f2c45f1SShri Abhyankar 4655f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This 4665f2c45f1SShri Abhyankar function is called during DMSetUp() */ 4675f2c45f1SShri Abhyankar #undef __FUNCT__ 4685f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkComponentSetUp" 4695f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 4705f2c45f1SShri Abhyankar { 4715f2c45f1SShri Abhyankar PetscErrorCode ierr; 4725f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4735f2c45f1SShri Abhyankar PetscInt arr_size; 4745f2c45f1SShri Abhyankar PetscInt p,offset,offsetp; 4755f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 4765f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 4775f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 4785f2c45f1SShri Abhyankar PetscInt ncomp, i; 4795f2c45f1SShri Abhyankar 4805f2c45f1SShri Abhyankar PetscFunctionBegin; 4815f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 4825f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 483*75b160a0SShri Abhyankar ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 4845f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 4855f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 4865f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 4875f2c45f1SShri Abhyankar /* Copy header */ 4885f2c45f1SShri Abhyankar header = &network->header[p]; 4895f2c45f1SShri Abhyankar ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType)); 4905f2c45f1SShri Abhyankar /* Copy data */ 4915f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 4925f2c45f1SShri Abhyankar ncomp = header->ndata; 4935f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 4945f2c45f1SShri Abhyankar offset = offsetp + network->dataheadersize + header->offset[i]; 4955f2c45f1SShri Abhyankar ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType)); 4965f2c45f1SShri Abhyankar } 4975f2c45f1SShri Abhyankar } 4985f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4995f2c45f1SShri Abhyankar } 5005f2c45f1SShri Abhyankar 5015f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 5025f2c45f1SShri Abhyankar #undef __FUNCT__ 5035f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkVariablesSetUp" 5045f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm) 5055f2c45f1SShri Abhyankar { 5065f2c45f1SShri Abhyankar PetscErrorCode ierr; 5075f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5085f2c45f1SShri Abhyankar 5095f2c45f1SShri Abhyankar PetscFunctionBegin; 5105f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 5115f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5125f2c45f1SShri Abhyankar } 5135f2c45f1SShri Abhyankar 5145f2c45f1SShri Abhyankar #undef __FUNCT__ 5155f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetComponentDataArray" 5165f2c45f1SShri Abhyankar /*@C 5175f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray - Returns the component data array 5185f2c45f1SShri Abhyankar 5195f2c45f1SShri Abhyankar Not Collective 5205f2c45f1SShri Abhyankar 5215f2c45f1SShri Abhyankar Input Parameters: 5225f2c45f1SShri Abhyankar . dm - The DMNetwork Object 5235f2c45f1SShri Abhyankar 5245f2c45f1SShri Abhyankar Output Parameters: 5255f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components 5265f2c45f1SShri Abhyankar 5275f2c45f1SShri Abhyankar Level: intermediate 5285f2c45f1SShri Abhyankar 5295f2c45f1SShri Abhyankar .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents 5305f2c45f1SShri Abhyankar @*/ 5315f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 5325f2c45f1SShri Abhyankar { 5335f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5345f2c45f1SShri Abhyankar 5355f2c45f1SShri Abhyankar PetscFunctionBegin; 5365f2c45f1SShri Abhyankar *componentdataarray = network->componentdataarray; 5375f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5385f2c45f1SShri Abhyankar } 5395f2c45f1SShri Abhyankar 5405f2c45f1SShri Abhyankar #undef __FUNCT__ 5415f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkDistribute" 5425f2c45f1SShri Abhyankar /*@ 5435f2c45f1SShri Abhyankar DMNetworkDistribute - Distributes the network and moves associated component data. 5445f2c45f1SShri Abhyankar 5455f2c45f1SShri Abhyankar Collective 5465f2c45f1SShri Abhyankar 5475f2c45f1SShri Abhyankar Input Parameter: 5485f2c45f1SShri Abhyankar + oldDM - the original DMNetwork object 5495f2c45f1SShri Abhyankar . partitioner - The partitioning package, or NULL for the default 5505f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default 5515f2c45f1SShri Abhyankar 5525f2c45f1SShri Abhyankar Output Parameter: 5535f2c45f1SShri Abhyankar . distDM - the distributed DMNetwork object 5545f2c45f1SShri Abhyankar 5555f2c45f1SShri Abhyankar Notes: 5565f2c45f1SShri Abhyankar This routine should be called only when using multiple processors. 5575f2c45f1SShri Abhyankar 5585f2c45f1SShri Abhyankar Distributes the network with a non-overlapping partitioning of the edges. 5595f2c45f1SShri Abhyankar 5605f2c45f1SShri Abhyankar Level: intermediate 5615f2c45f1SShri Abhyankar 5625f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 5635f2c45f1SShri Abhyankar @*/ 5645f2c45f1SShri Abhyankar PetscErrorCode DMNetworkDistribute(DM oldDM, const char partitioner[], PetscInt overlap,DM *distDM) 5655f2c45f1SShri Abhyankar { 5665f2c45f1SShri Abhyankar PetscErrorCode ierr; 5675f2c45f1SShri Abhyankar DM_Network *oldDMnetwork = (DM_Network*)oldDM->data; 5685f2c45f1SShri Abhyankar PetscSF pointsf; 5695f2c45f1SShri Abhyankar DM newDM; 5705f2c45f1SShri Abhyankar DM_Network *newDMnetwork; 5715f2c45f1SShri Abhyankar 5725f2c45f1SShri Abhyankar PetscFunctionBegin; 5735f2c45f1SShri Abhyankar ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr); 5745f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 5755f2c45f1SShri Abhyankar newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 5765f2c45f1SShri Abhyankar /* Distribute plex dm and dof section */ 5775f2c45f1SShri Abhyankar ierr = DMPlexDistribute(oldDMnetwork->plex,partitioner,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 5785f2c45f1SShri Abhyankar /* Distribute dof section */ 5795f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr); 5805f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 5815f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr); 5825f2c45f1SShri Abhyankar /* Distribute data and associated section */ 5835f2c45f1SShri Abhyankar ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPI_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 5845f2c45f1SShri Abhyankar /* Destroy point SF */ 5855f2c45f1SShri Abhyankar ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 5865f2c45f1SShri Abhyankar 5875f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 5885f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 5895f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 5905f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 5915f2c45f1SShri Abhyankar newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart; 5925f2c45f1SShri Abhyankar newDMnetwork->NNodes = oldDMnetwork->NNodes; 5935f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 5945f2c45f1SShri Abhyankar /* Set Dof section as the default section for dm */ 5955f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 5965f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 5975f2c45f1SShri Abhyankar 5985f2c45f1SShri Abhyankar *distDM = newDM; 5995f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6005f2c45f1SShri Abhyankar } 6015f2c45f1SShri Abhyankar 6025f2c45f1SShri Abhyankar #undef __FUNCT__ 6035f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetSupportingEdges" 6045f2c45f1SShri Abhyankar /*@C 6055f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 6065f2c45f1SShri Abhyankar 6075f2c45f1SShri Abhyankar Not Collective 6085f2c45f1SShri Abhyankar 6095f2c45f1SShri Abhyankar Input Parameters: 6105f2c45f1SShri Abhyankar + dm - The DMNetwork object 6115f2c45f1SShri Abhyankar - p - the vertex point 6125f2c45f1SShri Abhyankar 6135f2c45f1SShri Abhyankar Output Paramters: 6145f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 6155f2c45f1SShri Abhyankar - edges - List of edge points 6165f2c45f1SShri Abhyankar 6175f2c45f1SShri Abhyankar Level: intermediate 6185f2c45f1SShri Abhyankar 6195f2c45f1SShri Abhyankar Fortran Notes: 6205f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 6215f2c45f1SShri Abhyankar include petsc.h90 in your code. 6225f2c45f1SShri Abhyankar 6235f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes 6245f2c45f1SShri Abhyankar @*/ 6255f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 6265f2c45f1SShri Abhyankar { 6275f2c45f1SShri Abhyankar PetscErrorCode ierr; 6285f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 6295f2c45f1SShri Abhyankar 6305f2c45f1SShri Abhyankar PetscFunctionBegin; 6315f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 6325f2c45f1SShri Abhyankar ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 6335f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6345f2c45f1SShri Abhyankar } 6355f2c45f1SShri Abhyankar 6365f2c45f1SShri Abhyankar #undef __FUNCT__ 6375f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkGetConnectedNodes" 6385f2c45f1SShri Abhyankar /*@C 6395f2c45f1SShri Abhyankar DMNetworkGetConnectedNodes - Return the connected edges for this edge point 6405f2c45f1SShri Abhyankar 6415f2c45f1SShri Abhyankar Not Collective 6425f2c45f1SShri Abhyankar 6435f2c45f1SShri Abhyankar Input Parameters: 6445f2c45f1SShri Abhyankar + dm - The DMNetwork object 6455f2c45f1SShri Abhyankar - p - the edge point 6465f2c45f1SShri Abhyankar 6475f2c45f1SShri Abhyankar Output Paramters: 6485f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 6495f2c45f1SShri Abhyankar 6505f2c45f1SShri Abhyankar Level: intermediate 6515f2c45f1SShri Abhyankar 6525f2c45f1SShri Abhyankar Fortran Notes: 6535f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 6545f2c45f1SShri Abhyankar include petsc.h90 in your code. 6555f2c45f1SShri Abhyankar 6565f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 6575f2c45f1SShri Abhyankar @*/ 6585f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[]) 6595f2c45f1SShri Abhyankar { 6605f2c45f1SShri Abhyankar PetscErrorCode ierr; 6615f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 6625f2c45f1SShri Abhyankar 6635f2c45f1SShri Abhyankar PetscFunctionBegin; 6645f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 6655f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6665f2c45f1SShri Abhyankar } 6675f2c45f1SShri Abhyankar 6685f2c45f1SShri Abhyankar #undef __FUNCT__ 6695f2c45f1SShri Abhyankar #define __FUNCT__ "DMNetworkIsGhostVertex" 6705f2c45f1SShri Abhyankar /*@ 6715f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 6725f2c45f1SShri Abhyankar 6735f2c45f1SShri Abhyankar Not Collective 6745f2c45f1SShri Abhyankar 6755f2c45f1SShri Abhyankar Input Parameters: 6765f2c45f1SShri Abhyankar + dm - The DMNetwork object 6775f2c45f1SShri Abhyankar . p - the vertex point 6785f2c45f1SShri Abhyankar 6795f2c45f1SShri Abhyankar Output Parameter: 6805f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 6815f2c45f1SShri Abhyankar 6825f2c45f1SShri Abhyankar Level: intermediate 6835f2c45f1SShri Abhyankar 6845f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange 6855f2c45f1SShri Abhyankar @*/ 6865f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 6875f2c45f1SShri Abhyankar { 6885f2c45f1SShri Abhyankar PetscErrorCode ierr; 6895f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 6905f2c45f1SShri Abhyankar PetscInt offsetg; 6915f2c45f1SShri Abhyankar PetscSection sectiong; 6925f2c45f1SShri Abhyankar 6935f2c45f1SShri Abhyankar PetscFunctionBegin; 6945f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 6955f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 6965f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 6975f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 6985f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6995f2c45f1SShri Abhyankar } 7005f2c45f1SShri Abhyankar 7015f2c45f1SShri Abhyankar #undef __FUNCT__ 7025f2c45f1SShri Abhyankar #define __FUNCT__ "DMSetUp_Network" 7035f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 7045f2c45f1SShri Abhyankar { 7055f2c45f1SShri Abhyankar PetscErrorCode ierr; 7065f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 7075f2c45f1SShri Abhyankar 7085f2c45f1SShri Abhyankar PetscFunctionBegin; 7095f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 7105f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 7115f2c45f1SShri Abhyankar 7125f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 7135f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 7145f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7155f2c45f1SShri Abhyankar } 7165f2c45f1SShri Abhyankar 7175f2c45f1SShri Abhyankar #undef __FUNCT__ 7185f2c45f1SShri Abhyankar #define __FUNCT__ "DMCreateMatrix_Network" 7195f2c45f1SShri Abhyankar PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 7205f2c45f1SShri Abhyankar { 7215f2c45f1SShri Abhyankar PetscErrorCode ierr; 7225f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 7235f2c45f1SShri Abhyankar 7245f2c45f1SShri Abhyankar PetscFunctionBegin; 7255f2c45f1SShri Abhyankar ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 7265f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 7275f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7285f2c45f1SShri Abhyankar } 7295f2c45f1SShri Abhyankar 7305f2c45f1SShri Abhyankar #undef __FUNCT__ 7315f2c45f1SShri Abhyankar #define __FUNCT__ "DMDestroy_Network" 7325f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 7335f2c45f1SShri Abhyankar { 7345f2c45f1SShri Abhyankar PetscErrorCode ierr; 7355f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 7365f2c45f1SShri Abhyankar 7375f2c45f1SShri Abhyankar PetscFunctionBegin; 7385f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 7395f2c45f1SShri Abhyankar network->edges = NULL; 7405f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 7415f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 7425f2c45f1SShri Abhyankar /* ierr = PetscSectionDestroy(&network->GlobalDofSection);CHKERRQ(ierr); */ 7435f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 7445f2c45f1SShri Abhyankar ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 7455f2c45f1SShri Abhyankar ierr = PetscFree(network->header);CHKERRQ(ierr); 7465f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 7475f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7485f2c45f1SShri Abhyankar } 7495f2c45f1SShri Abhyankar 7505f2c45f1SShri Abhyankar #undef __FUNCT__ 7515f2c45f1SShri Abhyankar #define __FUNCT__ "DMView_Network" 7525f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 7535f2c45f1SShri Abhyankar { 7545f2c45f1SShri Abhyankar PetscErrorCode ierr; 7555f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 7565f2c45f1SShri Abhyankar 7575f2c45f1SShri Abhyankar PetscFunctionBegin; 7585f2c45f1SShri Abhyankar ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 7595f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7605f2c45f1SShri Abhyankar } 7615f2c45f1SShri Abhyankar 7625f2c45f1SShri Abhyankar #undef __FUNCT__ 7635f2c45f1SShri Abhyankar #define __FUNCT__ "DMGlobalToLocalBegin_Network" 7645f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 7655f2c45f1SShri Abhyankar { 7665f2c45f1SShri Abhyankar PetscErrorCode ierr; 7675f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 7685f2c45f1SShri Abhyankar 7695f2c45f1SShri Abhyankar PetscFunctionBegin; 7705f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 7715f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7725f2c45f1SShri Abhyankar } 7735f2c45f1SShri Abhyankar 7745f2c45f1SShri Abhyankar #undef __FUNCT__ 7755f2c45f1SShri Abhyankar #define __FUNCT__ "DMGlobalToLocalEnd_Network" 7765f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 7775f2c45f1SShri Abhyankar { 7785f2c45f1SShri Abhyankar PetscErrorCode ierr; 7795f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 7805f2c45f1SShri Abhyankar 7815f2c45f1SShri Abhyankar PetscFunctionBegin; 7825f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 7835f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7845f2c45f1SShri Abhyankar } 7855f2c45f1SShri Abhyankar 7865f2c45f1SShri Abhyankar #undef __FUNCT__ 7875f2c45f1SShri Abhyankar #define __FUNCT__ "DMLocalToGlobalBegin_Network" 7885f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 7895f2c45f1SShri Abhyankar { 7905f2c45f1SShri Abhyankar PetscErrorCode ierr; 7915f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 7925f2c45f1SShri Abhyankar 7935f2c45f1SShri Abhyankar PetscFunctionBegin; 7945f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 7955f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7965f2c45f1SShri Abhyankar } 7975f2c45f1SShri Abhyankar 7985f2c45f1SShri Abhyankar #undef __FUNCT__ 7995f2c45f1SShri Abhyankar #define __FUNCT__ "DMLocalToGlobalEnd_Network" 8005f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 8015f2c45f1SShri Abhyankar { 8025f2c45f1SShri Abhyankar PetscErrorCode ierr; 8035f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 8045f2c45f1SShri Abhyankar 8055f2c45f1SShri Abhyankar PetscFunctionBegin; 8065f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 8075f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8085f2c45f1SShri Abhyankar } 809