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 /*@ 65f2c45f1SShri Abhyankar DMNetworkSetSizes - Sets the local and global vertices and edges. 75f2c45f1SShri Abhyankar 85f2c45f1SShri Abhyankar Collective on DM 95f2c45f1SShri Abhyankar 105f2c45f1SShri Abhyankar Input Parameters: 115f2c45f1SShri Abhyankar + dm - the dm object 125f2c45f1SShri Abhyankar . nV - number of local vertices 135f2c45f1SShri Abhyankar . nE - number of local edges 145f2c45f1SShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) 155f2c45f1SShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE) 165f2c45f1SShri Abhyankar 175f2c45f1SShri Abhyankar Notes 185f2c45f1SShri Abhyankar If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 195f2c45f1SShri Abhyankar 205f2c45f1SShri Abhyankar You cannot change the sizes once they have been set 215f2c45f1SShri Abhyankar 221b266c99SBarry Smith Level: intermediate 231b266c99SBarry Smith 241b266c99SBarry Smith .seealso: DMNetworkCreate() 255f2c45f1SShri Abhyankar @*/ 265f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE) 275f2c45f1SShri Abhyankar { 285f2c45f1SShri Abhyankar PetscErrorCode ierr; 295f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 305f2c45f1SShri Abhyankar PetscInt a[2],b[2]; 315f2c45f1SShri Abhyankar 325f2c45f1SShri Abhyankar PetscFunctionBegin; 335f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 345f2c45f1SShri Abhyankar if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4); 355f2c45f1SShri Abhyankar if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5); 365f2c45f1SShri 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); 375f2c45f1SShri 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); 385f2c45f1SShri 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); 395f2c45f1SShri 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); 405f2c45f1SShri Abhyankar if (NE < 0 || NV < 0) { 415f2c45f1SShri Abhyankar a[0] = nV; a[1] = nE; 42b2566f29SBarry Smith ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 435f2c45f1SShri Abhyankar NV = b[0]; NE = b[1]; 445f2c45f1SShri Abhyankar } 455f2c45f1SShri Abhyankar network->nNodes = nV; 465f2c45f1SShri Abhyankar network->NNodes = NV; 475f2c45f1SShri Abhyankar network->nEdges = nE; 485f2c45f1SShri Abhyankar network->NEdges = NE; 495f2c45f1SShri Abhyankar PetscFunctionReturn(0); 505f2c45f1SShri Abhyankar } 515f2c45f1SShri Abhyankar 525f2c45f1SShri Abhyankar /*@ 535f2c45f1SShri Abhyankar DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 545f2c45f1SShri Abhyankar 555f2c45f1SShri Abhyankar Logically collective on DM 565f2c45f1SShri Abhyankar 575f2c45f1SShri Abhyankar Input Parameters: 585f2c45f1SShri Abhyankar . edges - list of edges 595f2c45f1SShri Abhyankar 605f2c45f1SShri Abhyankar Notes: 615f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 625f2c45f1SShri Abhyankar not be destroyed before the call to DMNetworkLayoutSetUp 635f2c45f1SShri Abhyankar 645f2c45f1SShri Abhyankar Level: intermediate 655f2c45f1SShri Abhyankar 665f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes 675f2c45f1SShri Abhyankar @*/ 685f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[]) 695f2c45f1SShri Abhyankar { 705f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 715f2c45f1SShri Abhyankar 725f2c45f1SShri Abhyankar PetscFunctionBegin; 735f2c45f1SShri Abhyankar network->edges = edgelist; 745f2c45f1SShri Abhyankar PetscFunctionReturn(0); 755f2c45f1SShri Abhyankar } 765f2c45f1SShri Abhyankar 775f2c45f1SShri Abhyankar /*@ 785f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 795f2c45f1SShri Abhyankar 805f2c45f1SShri Abhyankar Collective on DM 815f2c45f1SShri Abhyankar 825f2c45f1SShri Abhyankar Input Parameters 835f2c45f1SShri Abhyankar . DM - the dmnetwork object 845f2c45f1SShri Abhyankar 855f2c45f1SShri Abhyankar Notes: 865f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 875f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 885f2c45f1SShri Abhyankar 895f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 905f2c45f1SShri Abhyankar 915f2c45f1SShri Abhyankar Level: intermediate 925f2c45f1SShri Abhyankar 935f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 945f2c45f1SShri Abhyankar @*/ 955f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 965f2c45f1SShri Abhyankar { 975f2c45f1SShri Abhyankar PetscErrorCode ierr; 985f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 995f2c45f1SShri Abhyankar PetscInt dim = 1; /* One dimensional network */ 1005f2c45f1SShri Abhyankar PetscInt numCorners=2; 1015f2c45f1SShri Abhyankar PetscInt spacedim=2; 1025f2c45f1SShri Abhyankar double *vertexcoords=NULL; 1035f2c45f1SShri Abhyankar PetscInt i; 1045f2c45f1SShri Abhyankar PetscInt ndata; 1055f2c45f1SShri Abhyankar 1065f2c45f1SShri Abhyankar PetscFunctionBegin; 1075f2c45f1SShri Abhyankar if (network->nNodes) { 1089045477aSSatish Balay ierr = PetscCalloc1(numCorners*network->nNodes,&vertexcoords);CHKERRQ(ierr); 1095f2c45f1SShri Abhyankar } 1105f2c45f1SShri Abhyankar ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 1115f2c45f1SShri Abhyankar if (network->nNodes) { 1125f2c45f1SShri Abhyankar ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 1135f2c45f1SShri Abhyankar } 1145f2c45f1SShri Abhyankar ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 1155f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 1165f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 1175f2c45f1SShri Abhyankar 1185f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 1195f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 1205f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 1215f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 1225f2c45f1SShri Abhyankar 12324121865SAdrian Maldonado 12424121865SAdrian Maldonado 1255f2c45f1SShri Abhyankar network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 1266caa05f4SBarry Smith ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 1275f2c45f1SShri Abhyankar for (i = network->pStart; i < network->pEnd; i++) { 1285f2c45f1SShri Abhyankar network->header[i].ndata = 0; 1295f2c45f1SShri Abhyankar ndata = network->header[i].ndata; 1305f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 1315f2c45f1SShri Abhyankar network->header[i].offset[ndata] = 0; 1325f2c45f1SShri Abhyankar } 133854ce69bSBarry Smith ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 1345f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1355f2c45f1SShri Abhyankar } 1365f2c45f1SShri Abhyankar 137*94ef8ddeSSatish Balay /*@C 1385f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 1395f2c45f1SShri Abhyankar 1405f2c45f1SShri Abhyankar Logically collective on DM 1415f2c45f1SShri Abhyankar 1425f2c45f1SShri Abhyankar Input Parameters 1435f2c45f1SShri Abhyankar + dm - the network object 1445f2c45f1SShri Abhyankar . name - the component name 1455f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 1465f2c45f1SShri Abhyankar 1475f2c45f1SShri Abhyankar Output Parameters 1485f2c45f1SShri Abhyankar . key - an integer key that defines the component 1495f2c45f1SShri Abhyankar 1505f2c45f1SShri Abhyankar Notes 1515f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 1525f2c45f1SShri Abhyankar 1535f2c45f1SShri Abhyankar Level: intermediate 1545f2c45f1SShri Abhyankar 1555f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 1565f2c45f1SShri Abhyankar @*/ 1575f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 1585f2c45f1SShri Abhyankar { 1595f2c45f1SShri Abhyankar PetscErrorCode ierr; 1605f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1615f2c45f1SShri Abhyankar DMNetworkComponent *component=&network->component[network->ncomponent]; 1625f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 1635f2c45f1SShri Abhyankar PetscInt i; 1645f2c45f1SShri Abhyankar 1655f2c45f1SShri Abhyankar PetscFunctionBegin; 1665f2c45f1SShri Abhyankar 1675f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 1685f2c45f1SShri Abhyankar ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 1695f2c45f1SShri Abhyankar if (flg) { 1705f2c45f1SShri Abhyankar *key = i; 1715f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1725f2c45f1SShri Abhyankar } 1735f2c45f1SShri Abhyankar } 1745f2c45f1SShri Abhyankar 1755f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 1765f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 1775f2c45f1SShri Abhyankar *key = network->ncomponent; 1785f2c45f1SShri Abhyankar network->ncomponent++; 1795f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1805f2c45f1SShri Abhyankar } 1815f2c45f1SShri Abhyankar 1825f2c45f1SShri Abhyankar /*@ 1835f2c45f1SShri Abhyankar DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 1845f2c45f1SShri Abhyankar 1855f2c45f1SShri Abhyankar Not Collective 1865f2c45f1SShri Abhyankar 1875f2c45f1SShri Abhyankar Input Parameters: 1885f2c45f1SShri Abhyankar + dm - The DMNetwork object 1895f2c45f1SShri Abhyankar 1905f2c45f1SShri Abhyankar Output Paramters: 1915f2c45f1SShri Abhyankar + vStart - The first vertex point 1925f2c45f1SShri Abhyankar - vEnd - One beyond the last vertex point 1935f2c45f1SShri Abhyankar 1945f2c45f1SShri Abhyankar Level: intermediate 1955f2c45f1SShri Abhyankar 1965f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange 1975f2c45f1SShri Abhyankar @*/ 1985f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 1995f2c45f1SShri Abhyankar { 2005f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2015f2c45f1SShri Abhyankar 2025f2c45f1SShri Abhyankar PetscFunctionBegin; 2035f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 2045f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 2055f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2065f2c45f1SShri Abhyankar } 2075f2c45f1SShri Abhyankar 2085f2c45f1SShri Abhyankar /*@ 2095f2c45f1SShri Abhyankar DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 2105f2c45f1SShri Abhyankar 2115f2c45f1SShri Abhyankar Not Collective 2125f2c45f1SShri Abhyankar 2135f2c45f1SShri Abhyankar Input Parameters: 2145f2c45f1SShri Abhyankar + dm - The DMNetwork object 2155f2c45f1SShri Abhyankar 2165f2c45f1SShri Abhyankar Output Paramters: 2175f2c45f1SShri Abhyankar + eStart - The first edge point 2185f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 2195f2c45f1SShri Abhyankar 2205f2c45f1SShri Abhyankar Level: intermediate 2215f2c45f1SShri Abhyankar 2225f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange 2235f2c45f1SShri Abhyankar @*/ 2245f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 2255f2c45f1SShri Abhyankar { 2265f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2275f2c45f1SShri Abhyankar 2285f2c45f1SShri Abhyankar PetscFunctionBegin; 2295f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 2305f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 2315f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2325f2c45f1SShri Abhyankar } 2335f2c45f1SShri Abhyankar 2345f2c45f1SShri Abhyankar /*@ 2355f2c45f1SShri Abhyankar DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 2365f2c45f1SShri Abhyankar 2375f2c45f1SShri Abhyankar Not Collective 2385f2c45f1SShri Abhyankar 2395f2c45f1SShri Abhyankar Input Parameters: 2405f2c45f1SShri Abhyankar + dm - The DMNetwork object 2415f2c45f1SShri Abhyankar . p - vertex/edge point 2425f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component 2435f2c45f1SShri Abhyankar - compvalue - pointer to the data structure for the component 2445f2c45f1SShri Abhyankar 2455f2c45f1SShri Abhyankar Level: intermediate 2465f2c45f1SShri Abhyankar 2475f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 2485f2c45f1SShri Abhyankar @*/ 2495f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 2505f2c45f1SShri Abhyankar { 2515f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 25243a39a44SBarry Smith DMNetworkComponent *component = &network->component[componentkey]; 2535f2c45f1SShri Abhyankar DMNetworkComponentHeader header = &network->header[p]; 2545f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue = &network->cvalue[p]; 2555f2c45f1SShri Abhyankar PetscErrorCode ierr; 2565f2c45f1SShri Abhyankar 2575f2c45f1SShri Abhyankar PetscFunctionBegin; 258fa58f0a9SHong 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); 259fa58f0a9SHong Zhang 26043a39a44SBarry Smith header->size[header->ndata] = component->size; 26143a39a44SBarry Smith ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 2625f2c45f1SShri Abhyankar header->key[header->ndata] = componentkey; 2635f2c45f1SShri Abhyankar if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 2645f2c45f1SShri Abhyankar 2655f2c45f1SShri Abhyankar cvalue->data[header->ndata] = (void*)compvalue; 2665f2c45f1SShri Abhyankar header->ndata++; 2675f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2685f2c45f1SShri Abhyankar } 2695f2c45f1SShri Abhyankar 2705f2c45f1SShri Abhyankar /*@ 2715f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 2725f2c45f1SShri Abhyankar 2735f2c45f1SShri Abhyankar Not Collective 2745f2c45f1SShri Abhyankar 2755f2c45f1SShri Abhyankar Input Parameters: 2765f2c45f1SShri Abhyankar + dm - The DMNetwork object 2775f2c45f1SShri Abhyankar . p - vertex/edge point 2785f2c45f1SShri Abhyankar 2795f2c45f1SShri Abhyankar Output Parameters: 2805f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 2815f2c45f1SShri Abhyankar 2825f2c45f1SShri Abhyankar Level: intermediate 2835f2c45f1SShri Abhyankar 2845f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 2855f2c45f1SShri Abhyankar @*/ 2865f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 2875f2c45f1SShri Abhyankar { 2885f2c45f1SShri Abhyankar PetscErrorCode ierr; 2895f2c45f1SShri Abhyankar PetscInt offset; 2905f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2915f2c45f1SShri Abhyankar 2925f2c45f1SShri Abhyankar PetscFunctionBegin; 2935f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 2945f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 2955f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2965f2c45f1SShri Abhyankar } 2975f2c45f1SShri Abhyankar 2985f2c45f1SShri Abhyankar /*@ 2995f2c45f1SShri Abhyankar DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the 3005f2c45f1SShri Abhyankar component value from the component data array 3015f2c45f1SShri Abhyankar 3025f2c45f1SShri Abhyankar Not Collective 3035f2c45f1SShri Abhyankar 3045f2c45f1SShri Abhyankar Input Parameters: 3055f2c45f1SShri Abhyankar + dm - The DMNetwork object 3065f2c45f1SShri Abhyankar . p - vertex/edge point 3075f2c45f1SShri Abhyankar - compnum - component number 3085f2c45f1SShri Abhyankar 3095f2c45f1SShri Abhyankar Output Parameters: 3105f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component 3115f2c45f1SShri Abhyankar - offset - offset into the component data array associated with the vertex/edge point 3125f2c45f1SShri Abhyankar 3135f2c45f1SShri Abhyankar Notes: 3145f2c45f1SShri Abhyankar Typical usage: 3155f2c45f1SShri Abhyankar 3165f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray(dm, &arr); 3175f2c45f1SShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 3185f2c45f1SShri Abhyankar Loop over vertices or edges 3195f2c45f1SShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 3205f2c45f1SShri Abhyankar Loop over numcomps 3215f2c45f1SShri Abhyankar DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset); 3225f2c45f1SShri Abhyankar compdata = (UserCompDataType)(arr+offset); 3235f2c45f1SShri Abhyankar 3245f2c45f1SShri Abhyankar Level: intermediate 3255f2c45f1SShri Abhyankar 3265f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 3275f2c45f1SShri Abhyankar @*/ 3285f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 3295f2c45f1SShri Abhyankar { 3305f2c45f1SShri Abhyankar PetscErrorCode ierr; 3315f2c45f1SShri Abhyankar PetscInt offsetp; 3325f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 3335f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3345f2c45f1SShri Abhyankar 3355f2c45f1SShri Abhyankar PetscFunctionBegin; 3365f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 3375f2c45f1SShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 338d36f4e81SHong Zhang if (compkey) *compkey = header->key[compnum]; 339d36f4e81SHong Zhang if (offset) *offset = offsetp+network->dataheadersize+header->offset[compnum]; 3405f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3415f2c45f1SShri Abhyankar } 3425f2c45f1SShri Abhyankar 3435f2c45f1SShri Abhyankar /*@ 3445f2c45f1SShri Abhyankar DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 3455f2c45f1SShri Abhyankar 3465f2c45f1SShri Abhyankar Not Collective 3475f2c45f1SShri Abhyankar 3485f2c45f1SShri Abhyankar Input Parameters: 3495f2c45f1SShri Abhyankar + dm - The DMNetwork object 3505f2c45f1SShri Abhyankar - p - the edge/vertex point 3515f2c45f1SShri Abhyankar 3525f2c45f1SShri Abhyankar Output Parameters: 3535f2c45f1SShri Abhyankar . offset - the offset 3545f2c45f1SShri Abhyankar 3555f2c45f1SShri Abhyankar Level: intermediate 3565f2c45f1SShri Abhyankar 3575f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 3585f2c45f1SShri Abhyankar @*/ 3595f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 3605f2c45f1SShri Abhyankar { 3615f2c45f1SShri Abhyankar PetscErrorCode ierr; 3625f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3635f2c45f1SShri Abhyankar 3645f2c45f1SShri Abhyankar PetscFunctionBegin; 3655f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DofSection,p,offset);CHKERRQ(ierr); 3665f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3675f2c45f1SShri Abhyankar } 3685f2c45f1SShri Abhyankar 3695f2c45f1SShri Abhyankar /*@ 3705f2c45f1SShri Abhyankar DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 3715f2c45f1SShri Abhyankar 3725f2c45f1SShri Abhyankar Not Collective 3735f2c45f1SShri Abhyankar 3745f2c45f1SShri Abhyankar Input Parameters: 3755f2c45f1SShri Abhyankar + dm - The DMNetwork object 3765f2c45f1SShri Abhyankar - p - the edge/vertex point 3775f2c45f1SShri Abhyankar 3785f2c45f1SShri Abhyankar Output Parameters: 3795f2c45f1SShri Abhyankar . offsetg - the offset 3805f2c45f1SShri Abhyankar 3815f2c45f1SShri Abhyankar Level: intermediate 3825f2c45f1SShri Abhyankar 3835f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 3845f2c45f1SShri Abhyankar @*/ 3855f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 3865f2c45f1SShri Abhyankar { 3875f2c45f1SShri Abhyankar PetscErrorCode ierr; 3885f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3895f2c45f1SShri Abhyankar 3905f2c45f1SShri Abhyankar PetscFunctionBegin; 3915f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->GlobalDofSection,p,offsetg);CHKERRQ(ierr); 392dd6f46cdSHong Zhang if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost node */ 3935f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3945f2c45f1SShri Abhyankar } 3955f2c45f1SShri Abhyankar 39624121865SAdrian Maldonado /*@ 39724121865SAdrian Maldonado DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector. 39824121865SAdrian Maldonado 39924121865SAdrian Maldonado Not Collective 40024121865SAdrian Maldonado 40124121865SAdrian Maldonado Input Parameters: 40224121865SAdrian Maldonado + dm - The DMNetwork object 40324121865SAdrian Maldonado - p - the edge point 40424121865SAdrian Maldonado 40524121865SAdrian Maldonado Output Parameters: 40624121865SAdrian Maldonado . offset - the offset 40724121865SAdrian Maldonado 40824121865SAdrian Maldonado Level: intermediate 40924121865SAdrian Maldonado 41024121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 41124121865SAdrian Maldonado @*/ 41224121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset) 41324121865SAdrian Maldonado { 41424121865SAdrian Maldonado PetscErrorCode ierr; 41524121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 41624121865SAdrian Maldonado 41724121865SAdrian Maldonado PetscFunctionBegin; 41824121865SAdrian Maldonado 41924121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr); 42024121865SAdrian Maldonado PetscFunctionReturn(0); 42124121865SAdrian Maldonado } 42224121865SAdrian Maldonado 42324121865SAdrian Maldonado /*@ 42424121865SAdrian Maldonado DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector. 42524121865SAdrian Maldonado 42624121865SAdrian Maldonado Not Collective 42724121865SAdrian Maldonado 42824121865SAdrian Maldonado Input Parameters: 42924121865SAdrian Maldonado + dm - The DMNetwork object 43024121865SAdrian Maldonado - p - the vertex point 43124121865SAdrian Maldonado 43224121865SAdrian Maldonado Output Parameters: 43324121865SAdrian Maldonado . offset - the offset 43424121865SAdrian Maldonado 43524121865SAdrian Maldonado Level: intermediate 43624121865SAdrian Maldonado 43724121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 43824121865SAdrian Maldonado @*/ 43924121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset) 44024121865SAdrian Maldonado { 44124121865SAdrian Maldonado PetscErrorCode ierr; 44224121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 44324121865SAdrian Maldonado 44424121865SAdrian Maldonado PetscFunctionBegin; 44524121865SAdrian Maldonado 44624121865SAdrian Maldonado p -= network->vStart; 44724121865SAdrian Maldonado 44824121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr); 44924121865SAdrian Maldonado PetscFunctionReturn(0); 45024121865SAdrian Maldonado } 4515f2c45f1SShri Abhyankar /*@ 4525f2c45f1SShri Abhyankar DMNetworkAddNumVariables - Add number of variables associated with a given point. 4535f2c45f1SShri Abhyankar 4545f2c45f1SShri Abhyankar Not Collective 4555f2c45f1SShri Abhyankar 4565f2c45f1SShri Abhyankar Input Parameters: 4575f2c45f1SShri Abhyankar + dm - The DMNetworkObject 4585f2c45f1SShri Abhyankar . p - the vertex/edge point 4595f2c45f1SShri Abhyankar - nvar - number of additional variables 4605f2c45f1SShri Abhyankar 4615f2c45f1SShri Abhyankar Level: intermediate 4625f2c45f1SShri Abhyankar 4635f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables 4645f2c45f1SShri Abhyankar @*/ 4655f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 4665f2c45f1SShri Abhyankar { 4675f2c45f1SShri Abhyankar PetscErrorCode ierr; 4685f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4695f2c45f1SShri Abhyankar 4705f2c45f1SShri Abhyankar PetscFunctionBegin; 4715f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 4725f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4735f2c45f1SShri Abhyankar } 4745f2c45f1SShri Abhyankar 47527f51fceSHong Zhang /*@ 47627f51fceSHong Zhang DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point. 47727f51fceSHong Zhang 47827f51fceSHong Zhang Not Collective 47927f51fceSHong Zhang 48027f51fceSHong Zhang Input Parameters: 48127f51fceSHong Zhang + dm - The DMNetworkObject 48227f51fceSHong Zhang - p - the vertex/edge point 48327f51fceSHong Zhang 48427f51fceSHong Zhang Output Parameters: 48527f51fceSHong Zhang . nvar - number of variables 48627f51fceSHong Zhang 48727f51fceSHong Zhang Level: intermediate 48827f51fceSHong Zhang 48927f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables 49027f51fceSHong Zhang @*/ 49127f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar) 49227f51fceSHong Zhang { 49327f51fceSHong Zhang PetscErrorCode ierr; 49427f51fceSHong Zhang DM_Network *network = (DM_Network*)dm->data; 49527f51fceSHong Zhang 49627f51fceSHong Zhang PetscFunctionBegin; 49727f51fceSHong Zhang ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 49827f51fceSHong Zhang PetscFunctionReturn(0); 49927f51fceSHong Zhang } 50027f51fceSHong Zhang 5015f2c45f1SShri Abhyankar /*@ 5025f2c45f1SShri Abhyankar DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 5035f2c45f1SShri Abhyankar 5045f2c45f1SShri Abhyankar Not Collective 5055f2c45f1SShri Abhyankar 5065f2c45f1SShri Abhyankar Input Parameters: 5075f2c45f1SShri Abhyankar + dm - The DMNetworkObject 5085f2c45f1SShri Abhyankar . p - the vertex/edge point 5095f2c45f1SShri Abhyankar - nvar - number of variables 5105f2c45f1SShri Abhyankar 5115f2c45f1SShri Abhyankar Level: intermediate 5125f2c45f1SShri Abhyankar 5135f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables 5145f2c45f1SShri Abhyankar @*/ 5155f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 5165f2c45f1SShri Abhyankar { 5175f2c45f1SShri Abhyankar PetscErrorCode ierr; 5185f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5195f2c45f1SShri Abhyankar 5205f2c45f1SShri Abhyankar PetscFunctionBegin; 5215f2c45f1SShri Abhyankar ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 5225f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5235f2c45f1SShri Abhyankar } 5245f2c45f1SShri Abhyankar 5255f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This 5265f2c45f1SShri Abhyankar function is called during DMSetUp() */ 5275f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 5285f2c45f1SShri Abhyankar { 5295f2c45f1SShri Abhyankar PetscErrorCode ierr; 5305f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5315f2c45f1SShri Abhyankar PetscInt arr_size; 5325f2c45f1SShri Abhyankar PetscInt p,offset,offsetp; 5335f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 5345f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 5355f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 5365f2c45f1SShri Abhyankar PetscInt ncomp, i; 5375f2c45f1SShri Abhyankar 5385f2c45f1SShri Abhyankar PetscFunctionBegin; 5395f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 5405f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 54175b160a0SShri Abhyankar ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 5425f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 5435f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 5445f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 5455f2c45f1SShri Abhyankar /* Copy header */ 5465f2c45f1SShri Abhyankar header = &network->header[p]; 547302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 5485f2c45f1SShri Abhyankar /* Copy data */ 5495f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 5505f2c45f1SShri Abhyankar ncomp = header->ndata; 5515f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 5525f2c45f1SShri Abhyankar offset = offsetp + network->dataheadersize + header->offset[i]; 553302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 5545f2c45f1SShri Abhyankar } 5555f2c45f1SShri Abhyankar } 5565f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5575f2c45f1SShri Abhyankar } 5585f2c45f1SShri Abhyankar 5595f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 5605f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm) 5615f2c45f1SShri Abhyankar { 5625f2c45f1SShri Abhyankar PetscErrorCode ierr; 5635f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5645f2c45f1SShri Abhyankar 5655f2c45f1SShri Abhyankar PetscFunctionBegin; 5665f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 5675f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5685f2c45f1SShri Abhyankar } 5695f2c45f1SShri Abhyankar 5705f2c45f1SShri Abhyankar /*@C 5715f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray - Returns the component data array 5725f2c45f1SShri Abhyankar 5735f2c45f1SShri Abhyankar Not Collective 5745f2c45f1SShri Abhyankar 5755f2c45f1SShri Abhyankar Input Parameters: 5765f2c45f1SShri Abhyankar . dm - The DMNetwork Object 5775f2c45f1SShri Abhyankar 5785f2c45f1SShri Abhyankar Output Parameters: 5795f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components 5805f2c45f1SShri Abhyankar 5815f2c45f1SShri Abhyankar Level: intermediate 5825f2c45f1SShri Abhyankar 5835f2c45f1SShri Abhyankar .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents 5845f2c45f1SShri Abhyankar @*/ 5855f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 5865f2c45f1SShri Abhyankar { 5875f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5885f2c45f1SShri Abhyankar 5895f2c45f1SShri Abhyankar PetscFunctionBegin; 5905f2c45f1SShri Abhyankar *componentdataarray = network->componentdataarray; 5915f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5925f2c45f1SShri Abhyankar } 5935f2c45f1SShri Abhyankar 59424121865SAdrian Maldonado /* Get a subsection from a range of points */ 59524121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection) 59624121865SAdrian Maldonado { 59724121865SAdrian Maldonado PetscErrorCode ierr; 59824121865SAdrian Maldonado PetscInt i, nvar; 59924121865SAdrian Maldonado 60024121865SAdrian Maldonado PetscFunctionBegin; 60124121865SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr); 60224121865SAdrian Maldonado ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr); 60324121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 60424121865SAdrian Maldonado ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr); 60524121865SAdrian Maldonado ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr); 60624121865SAdrian Maldonado } 60724121865SAdrian Maldonado 60824121865SAdrian Maldonado ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr); 60924121865SAdrian Maldonado PetscFunctionReturn(0); 61024121865SAdrian Maldonado } 61124121865SAdrian Maldonado 61224121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */ 61324121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 61424121865SAdrian Maldonado { 61524121865SAdrian Maldonado PetscErrorCode ierr; 61624121865SAdrian Maldonado PetscInt i, *subpoints; 61724121865SAdrian Maldonado 61824121865SAdrian Maldonado PetscFunctionBegin; 61924121865SAdrian Maldonado /* Create index sets to map from "points" to "subpoints" */ 62024121865SAdrian Maldonado ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr); 62124121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 62224121865SAdrian Maldonado subpoints[i - pstart] = i; 62324121865SAdrian Maldonado } 624459726d8SSatish Balay ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr); 62524121865SAdrian Maldonado ierr = PetscFree(subpoints);CHKERRQ(ierr); 62624121865SAdrian Maldonado PetscFunctionReturn(0); 62724121865SAdrian Maldonado } 62824121865SAdrian Maldonado 62924121865SAdrian Maldonado /*@ 63024121865SAdrian Maldonado DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute. 63124121865SAdrian Maldonado 63224121865SAdrian Maldonado Collective 63324121865SAdrian Maldonado 63424121865SAdrian Maldonado Input Parameters: 63524121865SAdrian Maldonado . dm - The DMNetworkObject 63624121865SAdrian Maldonado 63724121865SAdrian Maldonado Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 63824121865SAdrian Maldonado 63924121865SAdrian Maldonado points = [0 1 2 3 4 5 6] 64024121865SAdrian Maldonado 64124121865SAdrian Maldonado where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]). 64224121865SAdrian Maldonado 64324121865SAdrian Maldonado With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 64424121865SAdrian Maldonado 64524121865SAdrian Maldonado Level: intermediate 64624121865SAdrian Maldonado 64724121865SAdrian Maldonado @*/ 64824121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 64924121865SAdrian Maldonado { 65024121865SAdrian Maldonado PetscErrorCode ierr; 65124121865SAdrian Maldonado MPI_Comm comm; 6529852e123SBarry Smith PetscMPIInt rank, size; 65324121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 65424121865SAdrian Maldonado 655eab1376dSHong Zhang PetscFunctionBegin; 65624121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 65724121865SAdrian Maldonado ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 6589852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 65924121865SAdrian Maldonado 66024121865SAdrian Maldonado /* Create maps for vertices and edges */ 66124121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 66224121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr); 66324121865SAdrian Maldonado 66424121865SAdrian Maldonado /* Create local sub-sections */ 66524121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr); 66624121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr); 66724121865SAdrian Maldonado 6689852e123SBarry Smith if (size > 1) { 66924121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 67024121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr); 67124121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr); 67224121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr); 67324121865SAdrian Maldonado } else { 67424121865SAdrian Maldonado /* create structures for vertex */ 67524121865SAdrian Maldonado ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr); 67624121865SAdrian Maldonado /* create structures for edge */ 67724121865SAdrian Maldonado ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr); 67824121865SAdrian Maldonado } 67924121865SAdrian Maldonado 68024121865SAdrian Maldonado 68124121865SAdrian Maldonado /* Add viewers */ 68224121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr); 68324121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr); 68424121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr); 68524121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr); 68624121865SAdrian Maldonado 68724121865SAdrian Maldonado PetscFunctionReturn(0); 68824121865SAdrian Maldonado } 6895f2c45f1SShri Abhyankar /*@ 6905f2c45f1SShri Abhyankar DMNetworkDistribute - Distributes the network and moves associated component data. 6915f2c45f1SShri Abhyankar 6925f2c45f1SShri Abhyankar Collective 6935f2c45f1SShri Abhyankar 6945f2c45f1SShri Abhyankar Input Parameter: 6955f2c45f1SShri Abhyankar + oldDM - the original DMNetwork object 6965f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default 6975f2c45f1SShri Abhyankar 6985f2c45f1SShri Abhyankar Output Parameter: 6995f2c45f1SShri Abhyankar . distDM - the distributed DMNetwork object 7005f2c45f1SShri Abhyankar 7015f2c45f1SShri Abhyankar Notes: 7025f2c45f1SShri Abhyankar This routine should be called only when using multiple processors. 7035f2c45f1SShri Abhyankar 7048b171c8eSHong Zhang Distributes the network with <overlap>-overlapping partitioning of the edges. 7055f2c45f1SShri Abhyankar 7065f2c45f1SShri Abhyankar Level: intermediate 7075f2c45f1SShri Abhyankar 7085f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 7095f2c45f1SShri Abhyankar @*/ 71080cf41d5SMatthew G. Knepley PetscErrorCode DMNetworkDistribute(DM oldDM, PetscInt overlap,DM *distDM) 7115f2c45f1SShri Abhyankar { 7125f2c45f1SShri Abhyankar PetscErrorCode ierr; 7135f2c45f1SShri Abhyankar DM_Network *oldDMnetwork = (DM_Network*)oldDM->data; 7145f2c45f1SShri Abhyankar PetscSF pointsf; 7155f2c45f1SShri Abhyankar DM newDM; 7165f2c45f1SShri Abhyankar DM_Network *newDMnetwork; 7175f2c45f1SShri Abhyankar 7185f2c45f1SShri Abhyankar PetscFunctionBegin; 7195f2c45f1SShri Abhyankar ierr = DMNetworkCreate(PetscObjectComm((PetscObject)oldDM),&newDM);CHKERRQ(ierr); 7205f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 7215f2c45f1SShri Abhyankar newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 7225f2c45f1SShri Abhyankar /* Distribute plex dm and dof section */ 72380cf41d5SMatthew G. Knepley ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 7245f2c45f1SShri Abhyankar /* Distribute dof section */ 7255f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DofSection);CHKERRQ(ierr); 7265f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 7275f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)oldDM),&newDMnetwork->DataSection);CHKERRQ(ierr); 7285f2c45f1SShri Abhyankar /* Distribute data and associated section */ 72931da1fc8SHong Zhang ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 73024121865SAdrian Maldonado 7315f2c45f1SShri Abhyankar 7325f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 7335f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 7345f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 7355f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 7365f2c45f1SShri Abhyankar newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart; 7375f2c45f1SShri Abhyankar newDMnetwork->NNodes = oldDMnetwork->NNodes; 7385f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 73924121865SAdrian Maldonado 7405f2c45f1SShri Abhyankar /* Set Dof section as the default section for dm */ 7415f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 7425f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 7435f2c45f1SShri Abhyankar 74424121865SAdrian Maldonado /* Destroy point SF */ 74524121865SAdrian Maldonado ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 74624121865SAdrian Maldonado 7475f2c45f1SShri Abhyankar *distDM = newDM; 7485f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7495f2c45f1SShri Abhyankar } 7505f2c45f1SShri Abhyankar 75124121865SAdrian Maldonado /*@C 75224121865SAdrian Maldonado PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering. 75324121865SAdrian Maldonado 75424121865SAdrian Maldonado Input Parameters: 75524121865SAdrian Maldonado + masterSF - the original SF structure 75624121865SAdrian Maldonado - map - a ISLocalToGlobal mapping that contains the subset of points 75724121865SAdrian Maldonado 75824121865SAdrian Maldonado Output Parameters: 75924121865SAdrian Maldonado . subSF - a subset of the masterSF for the desired subset. 76024121865SAdrian Maldonado */ 76124121865SAdrian Maldonado 76224121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) { 76324121865SAdrian Maldonado 76424121865SAdrian Maldonado PetscErrorCode ierr; 76524121865SAdrian Maldonado PetscInt nroots, nleaves, *ilocal_sub; 76624121865SAdrian Maldonado PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 76724121865SAdrian Maldonado PetscInt *local_points, *remote_points; 76824121865SAdrian Maldonado PetscSFNode *iremote_sub; 76924121865SAdrian Maldonado const PetscInt *ilocal; 77024121865SAdrian Maldonado const PetscSFNode *iremote; 77124121865SAdrian Maldonado 77224121865SAdrian Maldonado PetscFunctionBegin; 77324121865SAdrian Maldonado ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 77424121865SAdrian Maldonado 77524121865SAdrian Maldonado /* Look for leaves that pertain to the subset of points. Get the local ordering */ 77624121865SAdrian Maldonado ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr); 77724121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr); 77824121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 77924121865SAdrian Maldonado if (ilocal_map[i] != -1) nleaves_sub += 1; 78024121865SAdrian Maldonado } 78124121865SAdrian Maldonado /* Re-number ilocal with subset numbering. Need information from roots */ 78224121865SAdrian Maldonado ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr); 78324121865SAdrian Maldonado for (i = 0; i < nroots; i++) local_points[i] = i; 78424121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr); 78524121865SAdrian Maldonado ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 78624121865SAdrian Maldonado ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 78724121865SAdrian Maldonado /* Fill up graph using local (that is, local to the subset) numbering. */ 7884b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr); 7894b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr); 79024121865SAdrian Maldonado nleaves_sub = 0; 79124121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 79224121865SAdrian Maldonado if (ilocal_map[i] != -1) { 79324121865SAdrian Maldonado ilocal_sub[nleaves_sub] = ilocal_map[i]; 7944b70a8deSAdrian Maldonado iremote_sub[nleaves_sub].rank = iremote[i].rank; 79524121865SAdrian Maldonado iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 79624121865SAdrian Maldonado nleaves_sub += 1; 79724121865SAdrian Maldonado } 79824121865SAdrian Maldonado } 79924121865SAdrian Maldonado ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr); 80024121865SAdrian Maldonado ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr); 80124121865SAdrian Maldonado 80224121865SAdrian Maldonado /* Create new subSF */ 80324121865SAdrian Maldonado ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr); 80424121865SAdrian Maldonado ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr); 8054b70a8deSAdrian Maldonado ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr); 80624121865SAdrian Maldonado ierr = PetscFree(ilocal_map);CHKERRQ(ierr); 8074b70a8deSAdrian Maldonado ierr = PetscFree(iremote_sub);CHKERRQ(ierr); 80824121865SAdrian Maldonado 80924121865SAdrian Maldonado PetscFunctionReturn(0); 81024121865SAdrian Maldonado } 81124121865SAdrian Maldonado 81224121865SAdrian Maldonado 8135f2c45f1SShri Abhyankar /*@C 8145f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 8155f2c45f1SShri Abhyankar 8165f2c45f1SShri Abhyankar Not Collective 8175f2c45f1SShri Abhyankar 8185f2c45f1SShri Abhyankar Input Parameters: 8195f2c45f1SShri Abhyankar + dm - The DMNetwork object 8205f2c45f1SShri Abhyankar - p - the vertex point 8215f2c45f1SShri Abhyankar 8225f2c45f1SShri Abhyankar Output Paramters: 8235f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 8245f2c45f1SShri Abhyankar - edges - List of edge points 8255f2c45f1SShri Abhyankar 8265f2c45f1SShri Abhyankar Level: intermediate 8275f2c45f1SShri Abhyankar 8285f2c45f1SShri Abhyankar Fortran Notes: 8295f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 8305f2c45f1SShri Abhyankar include petsc.h90 in your code. 8315f2c45f1SShri Abhyankar 8325f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes 8335f2c45f1SShri Abhyankar @*/ 8345f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 8355f2c45f1SShri Abhyankar { 8365f2c45f1SShri Abhyankar PetscErrorCode ierr; 8375f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8385f2c45f1SShri Abhyankar 8395f2c45f1SShri Abhyankar PetscFunctionBegin; 8405f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 8415f2c45f1SShri Abhyankar ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 8425f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8435f2c45f1SShri Abhyankar } 8445f2c45f1SShri Abhyankar 8455f2c45f1SShri Abhyankar /*@C 846596e729fSHong Zhang DMNetworkGetConnectedNodes - Return the connected vertices for this edge point 8475f2c45f1SShri Abhyankar 8485f2c45f1SShri Abhyankar Not Collective 8495f2c45f1SShri Abhyankar 8505f2c45f1SShri Abhyankar Input Parameters: 8515f2c45f1SShri Abhyankar + dm - The DMNetwork object 8525f2c45f1SShri Abhyankar - p - the edge point 8535f2c45f1SShri Abhyankar 8545f2c45f1SShri Abhyankar Output Paramters: 8555f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 8565f2c45f1SShri Abhyankar 8575f2c45f1SShri Abhyankar Level: intermediate 8585f2c45f1SShri Abhyankar 8595f2c45f1SShri Abhyankar Fortran Notes: 8605f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 8615f2c45f1SShri Abhyankar include petsc.h90 in your code. 8625f2c45f1SShri Abhyankar 8635f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 8645f2c45f1SShri Abhyankar @*/ 8655f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[]) 8665f2c45f1SShri Abhyankar { 8675f2c45f1SShri Abhyankar PetscErrorCode ierr; 8685f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8695f2c45f1SShri Abhyankar 8705f2c45f1SShri Abhyankar PetscFunctionBegin; 8715f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 8725f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8735f2c45f1SShri Abhyankar } 8745f2c45f1SShri Abhyankar 8755f2c45f1SShri Abhyankar /*@ 8765f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 8775f2c45f1SShri Abhyankar 8785f2c45f1SShri Abhyankar Not Collective 8795f2c45f1SShri Abhyankar 8805f2c45f1SShri Abhyankar Input Parameters: 8815f2c45f1SShri Abhyankar + dm - The DMNetwork object 8825f2c45f1SShri Abhyankar . p - the vertex point 8835f2c45f1SShri Abhyankar 8845f2c45f1SShri Abhyankar Output Parameter: 8855f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 8865f2c45f1SShri Abhyankar 8875f2c45f1SShri Abhyankar Level: intermediate 8885f2c45f1SShri Abhyankar 8895f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange 8905f2c45f1SShri Abhyankar @*/ 8915f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 8925f2c45f1SShri Abhyankar { 8935f2c45f1SShri Abhyankar PetscErrorCode ierr; 8945f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8955f2c45f1SShri Abhyankar PetscInt offsetg; 8965f2c45f1SShri Abhyankar PetscSection sectiong; 8975f2c45f1SShri Abhyankar 8985f2c45f1SShri Abhyankar PetscFunctionBegin; 8995f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 9005f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 9015f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 9025f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 9035f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9045f2c45f1SShri Abhyankar } 9055f2c45f1SShri Abhyankar 9065f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 9075f2c45f1SShri Abhyankar { 9085f2c45f1SShri Abhyankar PetscErrorCode ierr; 9095f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 9105f2c45f1SShri Abhyankar 9115f2c45f1SShri Abhyankar PetscFunctionBegin; 9125f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 9135f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 9145f2c45f1SShri Abhyankar 9155f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 9165f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 9175f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9185f2c45f1SShri Abhyankar } 9195f2c45f1SShri Abhyankar 9201ad426b7SHong Zhang /*@ 92117df6e9eSHong Zhang DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 9221ad426b7SHong Zhang -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 9231ad426b7SHong Zhang 9241ad426b7SHong Zhang Collective 9251ad426b7SHong Zhang 9261ad426b7SHong Zhang Input Parameters: 92783b2e829SHong Zhang + dm - The DMNetwork object 92883b2e829SHong Zhang . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 92983b2e829SHong Zhang - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 9301ad426b7SHong Zhang 9311ad426b7SHong Zhang Level: intermediate 9321ad426b7SHong Zhang 9331ad426b7SHong Zhang @*/ 93483b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 9351ad426b7SHong Zhang { 9361ad426b7SHong Zhang DM_Network *network=(DM_Network*)dm->data; 9371ad426b7SHong Zhang 9381ad426b7SHong Zhang PetscFunctionBegin; 93983b2e829SHong Zhang network->userEdgeJacobian = eflg; 94083b2e829SHong Zhang network->userVertexJacobian = vflg; 9411ad426b7SHong Zhang PetscFunctionReturn(0); 9421ad426b7SHong Zhang } 9431ad426b7SHong Zhang 9441ad426b7SHong Zhang /*@ 94583b2e829SHong Zhang DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 94683b2e829SHong Zhang 94783b2e829SHong Zhang Not Collective 94883b2e829SHong Zhang 94983b2e829SHong Zhang Input Parameters: 95083b2e829SHong Zhang + dm - The DMNetwork object 95183b2e829SHong Zhang . p - the edge point 9523e97b6e8SHong Zhang - J - array (size = 3) of Jacobian submatrices for this edge point: 9533e97b6e8SHong Zhang J[0]: this edge 9543e97b6e8SHong Zhang J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes() 95583b2e829SHong Zhang 95683b2e829SHong Zhang Level: intermediate 95783b2e829SHong Zhang 95883b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix 95983b2e829SHong Zhang @*/ 96083b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 96183b2e829SHong Zhang { 96283b2e829SHong Zhang PetscErrorCode ierr; 96383b2e829SHong Zhang DM_Network *network=(DM_Network*)dm->data; 96483b2e829SHong Zhang 96583b2e829SHong Zhang PetscFunctionBegin; 966883e35e8SHong Zhang if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 96783b2e829SHong Zhang if (!network->Je) { 968883e35e8SHong Zhang ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 96983b2e829SHong Zhang } 970883e35e8SHong Zhang network->Je[3*p] = J[0]; 971883e35e8SHong Zhang network->Je[3*p+1] = J[1]; 972883e35e8SHong Zhang network->Je[3*p+2] = J[2]; 97383b2e829SHong Zhang PetscFunctionReturn(0); 97483b2e829SHong Zhang } 97583b2e829SHong Zhang 97683b2e829SHong Zhang /*@ 97776ddfea5SHong Zhang DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 9781ad426b7SHong Zhang 9791ad426b7SHong Zhang Not Collective 9801ad426b7SHong Zhang 9811ad426b7SHong Zhang Input Parameters: 9821ad426b7SHong Zhang + dm - The DMNetwork object 9831ad426b7SHong Zhang . p - the vertex point 9843e97b6e8SHong Zhang - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 9853e97b6e8SHong Zhang J[0]: this vertex 9863e97b6e8SHong Zhang J[1+2*i]: i-th supporting edge 9873e97b6e8SHong Zhang J[1+2*i+1]: i-th connected vertex 9881ad426b7SHong Zhang 9891ad426b7SHong Zhang Level: intermediate 9901ad426b7SHong Zhang 99183b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix 9921ad426b7SHong Zhang @*/ 993883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 9945f2c45f1SShri Abhyankar { 9955f2c45f1SShri Abhyankar PetscErrorCode ierr; 9965f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 997f4431b8cSHong Zhang PetscInt i,*vptr,nedges,vStart,vEnd; 998883e35e8SHong Zhang const PetscInt *edges; 9995f2c45f1SShri Abhyankar 10005f2c45f1SShri Abhyankar PetscFunctionBegin; 1001883e35e8SHong Zhang if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 1002883e35e8SHong Zhang 1003883e35e8SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1004f4431b8cSHong Zhang 100583b2e829SHong Zhang if (!network->Jv) { 1006f4431b8cSHong Zhang PetscInt nNodes = network->nNodes,nedges_total; 1007883e35e8SHong Zhang /* count nvertex_total */ 10083e97b6e8SHong Zhang nedges_total = 0; 1009883e35e8SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1010f4431b8cSHong Zhang ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr); 1011f4431b8cSHong Zhang 1012883e35e8SHong Zhang vptr[0] = 0; 1013f4431b8cSHong Zhang for (i=0; i<nNodes; i++) { 1014f4431b8cSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 1015883e35e8SHong Zhang nedges_total += nedges; 1016f4431b8cSHong Zhang vptr[i+1] = vptr[i] + 2*nedges + 1; 10171ad426b7SHong Zhang } 10183e97b6e8SHong Zhang 1019f4431b8cSHong Zhang ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr); 1020883e35e8SHong Zhang network->Jvptr = vptr; 1021883e35e8SHong Zhang } 1022883e35e8SHong Zhang 1023883e35e8SHong Zhang vptr = network->Jvptr; 10243e97b6e8SHong Zhang network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 10253e97b6e8SHong Zhang 10263e97b6e8SHong Zhang /* Set Jacobian for each supporting edge and connected vertex */ 1027883e35e8SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 1028883e35e8SHong Zhang for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 1029883e35e8SHong Zhang PetscFunctionReturn(0); 1030883e35e8SHong Zhang } 1031883e35e8SHong Zhang 1032e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 10335cf7da58SHong Zhang { 10345cf7da58SHong Zhang PetscErrorCode ierr; 10355cf7da58SHong Zhang PetscInt j; 10365cf7da58SHong Zhang PetscScalar val=(PetscScalar)ncols; 10375cf7da58SHong Zhang 10385cf7da58SHong Zhang PetscFunctionBegin; 10395cf7da58SHong Zhang if (!ghost) { 10405cf7da58SHong Zhang for (j=0; j<nrows; j++) { 10415cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 10425cf7da58SHong Zhang } 10435cf7da58SHong Zhang } else { 10445cf7da58SHong Zhang for (j=0; j<nrows; j++) { 10455cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 10465cf7da58SHong Zhang } 10475cf7da58SHong Zhang } 10485cf7da58SHong Zhang PetscFunctionReturn(0); 10495cf7da58SHong Zhang } 10505cf7da58SHong Zhang 1051e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 10525cf7da58SHong Zhang { 10535cf7da58SHong Zhang PetscErrorCode ierr; 10545cf7da58SHong Zhang PetscInt j,ncols_u; 10555cf7da58SHong Zhang PetscScalar val; 10565cf7da58SHong Zhang 10575cf7da58SHong Zhang PetscFunctionBegin; 10585cf7da58SHong Zhang if (!ghost) { 10595cf7da58SHong Zhang for (j=0; j<nrows; j++) { 10605cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 10615cf7da58SHong Zhang val = (PetscScalar)ncols_u; 10625cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 10635cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 10645cf7da58SHong Zhang } 10655cf7da58SHong Zhang } else { 10665cf7da58SHong Zhang for (j=0; j<nrows; j++) { 10675cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 10685cf7da58SHong Zhang val = (PetscScalar)ncols_u; 10695cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 10705cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 10715cf7da58SHong Zhang } 10725cf7da58SHong Zhang } 10735cf7da58SHong Zhang PetscFunctionReturn(0); 10745cf7da58SHong Zhang } 10755cf7da58SHong Zhang 1076e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 10775cf7da58SHong Zhang { 10785cf7da58SHong Zhang PetscErrorCode ierr; 10795cf7da58SHong Zhang 10805cf7da58SHong Zhang PetscFunctionBegin; 10815cf7da58SHong Zhang if (Ju) { 10825cf7da58SHong Zhang ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 10835cf7da58SHong Zhang } else { 10845cf7da58SHong Zhang ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 10855cf7da58SHong Zhang } 10865cf7da58SHong Zhang PetscFunctionReturn(0); 10875cf7da58SHong Zhang } 10885cf7da58SHong Zhang 1089e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1090883e35e8SHong Zhang { 1091883e35e8SHong Zhang PetscErrorCode ierr; 1092883e35e8SHong Zhang PetscInt j,*cols; 1093883e35e8SHong Zhang PetscScalar *zeros; 1094883e35e8SHong Zhang 1095883e35e8SHong Zhang PetscFunctionBegin; 1096883e35e8SHong Zhang ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 1097883e35e8SHong Zhang for (j=0; j<ncols; j++) cols[j] = j+ cstart; 1098883e35e8SHong Zhang ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 1099883e35e8SHong Zhang ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 11001ad426b7SHong Zhang PetscFunctionReturn(0); 11011ad426b7SHong Zhang } 1102a4e85ca8SHong Zhang 1103e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 11043e97b6e8SHong Zhang { 11053e97b6e8SHong Zhang PetscErrorCode ierr; 11063e97b6e8SHong Zhang PetscInt j,M,N,row,col,ncols_u; 11073e97b6e8SHong Zhang const PetscInt *cols; 11083e97b6e8SHong Zhang PetscScalar zero=0.0; 11093e97b6e8SHong Zhang 11103e97b6e8SHong Zhang PetscFunctionBegin; 11113e97b6e8SHong Zhang ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 11123e97b6e8SHong 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); 11133e97b6e8SHong Zhang 11143e97b6e8SHong Zhang for (row=0; row<nrows; row++) { 11153e97b6e8SHong Zhang ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 11163e97b6e8SHong Zhang for (j=0; j<ncols_u; j++) { 11173e97b6e8SHong Zhang col = cols[j] + cstart; 11183e97b6e8SHong Zhang ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 11193e97b6e8SHong Zhang } 11203e97b6e8SHong Zhang ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 11213e97b6e8SHong Zhang } 11223e97b6e8SHong Zhang PetscFunctionReturn(0); 11233e97b6e8SHong Zhang } 11241ad426b7SHong Zhang 1125e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1126a4e85ca8SHong Zhang { 1127a4e85ca8SHong Zhang PetscErrorCode ierr; 1128f4431b8cSHong Zhang 1129a4e85ca8SHong Zhang PetscFunctionBegin; 1130a4e85ca8SHong Zhang if (Ju) { 1131a4e85ca8SHong Zhang ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1132a4e85ca8SHong Zhang } else { 1133a4e85ca8SHong Zhang ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1134a4e85ca8SHong Zhang } 1135a4e85ca8SHong Zhang PetscFunctionReturn(0); 1136a4e85ca8SHong Zhang } 1137a4e85ca8SHong Zhang 113824121865SAdrian Maldonado 113924121865SAdrian 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. 114024121865SAdrian Maldonado */ 114124121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 114224121865SAdrian Maldonado { 114324121865SAdrian Maldonado PetscErrorCode ierr; 114424121865SAdrian Maldonado PetscInt i, size, dof; 114524121865SAdrian Maldonado PetscInt *glob2loc; 114624121865SAdrian Maldonado 114724121865SAdrian Maldonado PetscFunctionBegin; 114824121865SAdrian Maldonado ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr); 114924121865SAdrian Maldonado ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr); 115024121865SAdrian Maldonado 115124121865SAdrian Maldonado for (i = 0; i < size; i++) { 115224121865SAdrian Maldonado ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr); 115324121865SAdrian Maldonado dof = (dof >= 0) ? dof : -(dof + 1); 115424121865SAdrian Maldonado glob2loc[i] = dof; 115524121865SAdrian Maldonado } 115624121865SAdrian Maldonado 115724121865SAdrian Maldonado ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 115824121865SAdrian Maldonado #if 0 115924121865SAdrian Maldonado ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 116024121865SAdrian Maldonado #endif 116124121865SAdrian Maldonado PetscFunctionReturn(0); 116224121865SAdrian Maldonado } 116324121865SAdrian Maldonado 11641ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 11651ad426b7SHong Zhang { 11661ad426b7SHong Zhang PetscErrorCode ierr; 116724121865SAdrian Maldonado PetscMPIInt rank, size; 11681ad426b7SHong Zhang DM_Network *network = (DM_Network*) dm->data; 1169a4e85ca8SHong Zhang PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 1170840c2264SHong Zhang PetscInt cstart,ncols,j,e,v; 117124121865SAdrian Maldonado PetscBool ghost,ghost_vc,ghost2,isNest; 1172a4e85ca8SHong Zhang Mat Juser; 1173bfbc38dcSHong Zhang PetscSection sectionGlobal; 1174447d78afSSatish Balay PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 1175a4e85ca8SHong Zhang const PetscInt *edges,*cone; 11765cf7da58SHong Zhang MPI_Comm comm; 117724121865SAdrian Maldonado MatType mtype; 11785cf7da58SHong Zhang Vec vd_nz,vo_nz; 11795cf7da58SHong Zhang PetscInt *dnnz,*onnz; 11805cf7da58SHong Zhang PetscScalar *vdnz,*vonz; 11811ad426b7SHong Zhang 11821ad426b7SHong Zhang PetscFunctionBegin; 118324121865SAdrian Maldonado mtype = dm->mattype; 118424121865SAdrian Maldonado ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr); 118524121865SAdrian Maldonado 118624121865SAdrian Maldonado if (isNest) { 11870731d606SHong Zhang /* ierr = DMCreateMatrix_Network_Nest(); */ 118824121865SAdrian Maldonado PetscInt eDof, vDof; 118924121865SAdrian Maldonado Mat j11, j12, j21, j22, bA[2][2]; 119024121865SAdrian Maldonado ISLocalToGlobalMapping eISMap, vISMap; 119124121865SAdrian Maldonado 119224121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 119324121865SAdrian Maldonado ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 119424121865SAdrian Maldonado ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 119524121865SAdrian Maldonado 119624121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr); 119724121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr); 119824121865SAdrian Maldonado 119924121865SAdrian Maldonado ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr); 120024121865SAdrian Maldonado ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 120124121865SAdrian Maldonado ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr); 120224121865SAdrian Maldonado 120324121865SAdrian Maldonado ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr); 120424121865SAdrian Maldonado ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr); 120524121865SAdrian Maldonado ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr); 120624121865SAdrian Maldonado 120724121865SAdrian Maldonado ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr); 120824121865SAdrian Maldonado ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 120924121865SAdrian Maldonado ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr); 121024121865SAdrian Maldonado 121124121865SAdrian Maldonado ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr); 121224121865SAdrian Maldonado ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 121324121865SAdrian Maldonado ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr); 121424121865SAdrian Maldonado 121524121865SAdrian Maldonado bA[0][0] = j11; 121624121865SAdrian Maldonado bA[0][1] = j12; 121724121865SAdrian Maldonado bA[1][0] = j21; 121824121865SAdrian Maldonado bA[1][1] = j22; 121924121865SAdrian Maldonado 122024121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr); 122124121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr); 122224121865SAdrian Maldonado 122324121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr); 122424121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr); 122524121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr); 122624121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr); 122724121865SAdrian Maldonado 122824121865SAdrian Maldonado ierr = MatSetUp(j11);CHKERRQ(ierr); 122924121865SAdrian Maldonado ierr = MatSetUp(j12);CHKERRQ(ierr); 123024121865SAdrian Maldonado ierr = MatSetUp(j21);CHKERRQ(ierr); 123124121865SAdrian Maldonado ierr = MatSetUp(j22);CHKERRQ(ierr); 123224121865SAdrian Maldonado 123324121865SAdrian Maldonado ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr); 123424121865SAdrian Maldonado ierr = MatSetUp(*J);CHKERRQ(ierr); 123524121865SAdrian Maldonado ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr); 123624121865SAdrian Maldonado ierr = MatDestroy(&j11);CHKERRQ(ierr); 123724121865SAdrian Maldonado ierr = MatDestroy(&j12);CHKERRQ(ierr); 123824121865SAdrian Maldonado ierr = MatDestroy(&j21);CHKERRQ(ierr); 123924121865SAdrian Maldonado ierr = MatDestroy(&j22);CHKERRQ(ierr); 124024121865SAdrian Maldonado 124124121865SAdrian Maldonado ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 124224121865SAdrian Maldonado ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 124324121865SAdrian Maldonado 124424121865SAdrian Maldonado /* Free structures */ 124524121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr); 124624121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr); 124724121865SAdrian Maldonado 124824121865SAdrian Maldonado PetscFunctionReturn(0); 124924121865SAdrian Maldonado } else if (!network->userEdgeJacobian && !network->userVertexJacobian) { 1250a4e85ca8SHong Zhang /* user does not provide Jacobian blocks */ 1251bfbc38dcSHong Zhang ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 1252bfbc38dcSHong Zhang ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 12531ad426b7SHong Zhang PetscFunctionReturn(0); 12541ad426b7SHong Zhang } 12551ad426b7SHong Zhang 1256bfbc38dcSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 12572a945128SHong Zhang ierr = DMGetDefaultGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 1258bfbc38dcSHong Zhang ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 1259bfbc38dcSHong Zhang ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 12602a945128SHong Zhang 12612a945128SHong Zhang ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 12622a945128SHong Zhang ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 126389898e50SHong Zhang 126489898e50SHong Zhang /* (1) Set matrix preallocation */ 126589898e50SHong Zhang /*------------------------------*/ 1266840c2264SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1267840c2264SHong Zhang ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr); 1268840c2264SHong Zhang ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr); 1269840c2264SHong Zhang ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr); 1270840c2264SHong Zhang ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr); 1271840c2264SHong Zhang ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr); 1272840c2264SHong Zhang 127389898e50SHong Zhang /* Set preallocation for edges */ 127489898e50SHong Zhang /*-----------------------------*/ 1275840c2264SHong Zhang ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 1276840c2264SHong Zhang 1277bdcb62a2SHong Zhang ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr); 1278840c2264SHong Zhang for (e=eStart; e<eEnd; e++) { 1279840c2264SHong Zhang /* Get row indices */ 1280840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 1281840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 1282840c2264SHong Zhang if (nrows) { 1283840c2264SHong Zhang if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je"); 1284840c2264SHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 1285840c2264SHong Zhang 12865cf7da58SHong Zhang /* Set preallocation for conntected vertices */ 1287840c2264SHong Zhang ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr); 1288840c2264SHong Zhang for (v=0; v<2; v++) { 1289840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 1290840c2264SHong Zhang 1291840c2264SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 1292840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 12935cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1294840c2264SHong Zhang } 1295840c2264SHong Zhang 129689898e50SHong Zhang /* Set preallocation for edge self */ 1297840c2264SHong Zhang cstart = rstart; 1298840c2264SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 12995cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1300840c2264SHong Zhang } 1301840c2264SHong Zhang } 1302840c2264SHong Zhang 130389898e50SHong Zhang /* Set preallocation for vertices */ 130489898e50SHong Zhang /*--------------------------------*/ 1305840c2264SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1306840c2264SHong Zhang if (vEnd - vStart) { 1307840c2264SHong Zhang if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv"); 1308840c2264SHong Zhang vptr = network->Jvptr; 1309840c2264SHong Zhang if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr"); 1310840c2264SHong Zhang } 1311840c2264SHong Zhang 1312840c2264SHong Zhang for (v=vStart; v<vEnd; v++) { 1313840c2264SHong Zhang /* Get row indices */ 1314840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1315840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 1316840c2264SHong Zhang if (!nrows) continue; 1317840c2264SHong Zhang 1318bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1319bdcb62a2SHong Zhang if (ghost) { 1320bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1321bdcb62a2SHong Zhang } else { 1322bdcb62a2SHong Zhang rows_v = rows; 1323bdcb62a2SHong Zhang } 1324bdcb62a2SHong Zhang 1325bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1326840c2264SHong Zhang 1327840c2264SHong Zhang /* Get supporting edges and connected vertices */ 1328840c2264SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1329840c2264SHong Zhang 1330840c2264SHong Zhang for (e=0; e<nedges; e++) { 1331840c2264SHong Zhang /* Supporting edges */ 1332840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1333840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1334840c2264SHong Zhang 1335840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1336bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1337840c2264SHong Zhang 1338840c2264SHong Zhang /* Connected vertices */ 1339840c2264SHong Zhang ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr); 1340840c2264SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 1341840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr); 1342840c2264SHong Zhang 1343840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1344840c2264SHong Zhang 1345840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1346e102a522SHong Zhang if (ghost_vc||ghost) { 1347e102a522SHong Zhang ghost2 = PETSC_TRUE; 1348e102a522SHong Zhang } else { 1349e102a522SHong Zhang ghost2 = PETSC_FALSE; 1350e102a522SHong Zhang } 1351e102a522SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr); 1352840c2264SHong Zhang } 1353840c2264SHong Zhang 135489898e50SHong Zhang /* Set preallocation for vertex self */ 1355840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1356840c2264SHong Zhang if (!ghost) { 1357840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1358840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1359bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1360840c2264SHong Zhang } 1361bdcb62a2SHong Zhang if (ghost) { 1362bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1363bdcb62a2SHong Zhang } 1364840c2264SHong Zhang } 1365840c2264SHong Zhang 1366840c2264SHong Zhang ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr); 1367840c2264SHong Zhang ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr); 13685cf7da58SHong Zhang 13695cf7da58SHong Zhang ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr); 13705cf7da58SHong Zhang 13715cf7da58SHong Zhang ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr); 1372840c2264SHong Zhang ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr); 1373840c2264SHong Zhang 1374840c2264SHong Zhang ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr); 1375840c2264SHong Zhang ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr); 1376840c2264SHong Zhang for (j=0; j<localSize; j++) { 1377e102a522SHong Zhang dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 1378e102a522SHong Zhang onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 1379840c2264SHong Zhang } 1380840c2264SHong Zhang ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr); 1381840c2264SHong Zhang ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr); 1382840c2264SHong Zhang ierr = VecDestroy(&vd_nz);CHKERRQ(ierr); 1383840c2264SHong Zhang ierr = VecDestroy(&vo_nz);CHKERRQ(ierr); 1384840c2264SHong Zhang 13855cf7da58SHong Zhang ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr); 13865cf7da58SHong Zhang ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr); 13875cf7da58SHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 13885cf7da58SHong Zhang 13895cf7da58SHong Zhang ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr); 13905cf7da58SHong Zhang 139189898e50SHong Zhang /* (2) Set matrix entries for edges */ 139289898e50SHong Zhang /*----------------------------------*/ 13931ad426b7SHong Zhang for (e=eStart; e<eEnd; e++) { 1394bfbc38dcSHong Zhang /* Get row indices */ 13951ad426b7SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 139617df6e9eSHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 13974b976069SHong Zhang if (nrows) { 13984b976069SHong Zhang if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je"); 1399bdcb62a2SHong Zhang 140017df6e9eSHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 14011ad426b7SHong Zhang 1402bfbc38dcSHong Zhang /* Set matrix entries for conntected vertices */ 140317df6e9eSHong Zhang ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr); 1404bfbc38dcSHong Zhang for (v=0; v<2; v++) { 1405bfbc38dcSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 1406883e35e8SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 14073e97b6e8SHong Zhang 1408a4e85ca8SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 1409a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1410bfbc38dcSHong Zhang } 141117df6e9eSHong Zhang 1412bfbc38dcSHong Zhang /* Set matrix entries for edge self */ 14133e97b6e8SHong Zhang cstart = rstart; 1414a4e85ca8SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 1415a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 14161ad426b7SHong Zhang } 14174b976069SHong Zhang } 14181ad426b7SHong Zhang 1419bfbc38dcSHong Zhang /* Set matrix entries for vertices */ 142083b2e829SHong Zhang /*---------------------------------*/ 14211ad426b7SHong Zhang for (v=vStart; v<vEnd; v++) { 1422bfbc38dcSHong Zhang /* Get row indices */ 1423596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1424596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 14254b976069SHong Zhang if (!nrows) continue; 1426596e729fSHong Zhang 1427bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1428bdcb62a2SHong Zhang if (ghost) { 1429bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1430bdcb62a2SHong Zhang } else { 1431bdcb62a2SHong Zhang rows_v = rows; 1432bdcb62a2SHong Zhang } 1433bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1434596e729fSHong Zhang 1435bfbc38dcSHong Zhang /* Get supporting edges and connected vertices */ 1436596e729fSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1437596e729fSHong Zhang 1438596e729fSHong Zhang for (e=0; e<nedges; e++) { 1439bfbc38dcSHong Zhang /* Supporting edges */ 1440596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1441596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1442596e729fSHong Zhang 1443a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1444bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1445596e729fSHong Zhang 1446bfbc38dcSHong Zhang /* Connected vertices */ 144744aca652SHong Zhang ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr); 14482a945128SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 14492a945128SHong Zhang 145044aca652SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 145144aca652SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1452a4e85ca8SHong Zhang 1453a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1454bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1455596e729fSHong Zhang } 1456596e729fSHong Zhang 1457bfbc38dcSHong Zhang /* Set matrix entries for vertex self */ 14581ad426b7SHong Zhang if (!ghost) { 1459596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1460a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1461bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr); 1462bdcb62a2SHong Zhang } 1463bdcb62a2SHong Zhang if (ghost) { 1464bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1465bdcb62a2SHong Zhang } 14661ad426b7SHong Zhang } 1467a4e85ca8SHong Zhang ierr = PetscFree(rows);CHKERRQ(ierr); 1468bdcb62a2SHong Zhang 14691ad426b7SHong Zhang ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14701ad426b7SHong Zhang ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1471dd6f46cdSHong Zhang 14725f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 14735f2c45f1SShri Abhyankar PetscFunctionReturn(0); 14745f2c45f1SShri Abhyankar } 14755f2c45f1SShri Abhyankar 14765f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 14775f2c45f1SShri Abhyankar { 14785f2c45f1SShri Abhyankar PetscErrorCode ierr; 14795f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 14805f2c45f1SShri Abhyankar 14815f2c45f1SShri Abhyankar PetscFunctionBegin; 14828415c774SShri Abhyankar if (--network->refct > 0) PetscFunctionReturn(0); 148383b2e829SHong Zhang if (network->Je) { 148483b2e829SHong Zhang ierr = PetscFree(network->Je);CHKERRQ(ierr); 148583b2e829SHong Zhang } 148683b2e829SHong Zhang if (network->Jv) { 1487883e35e8SHong Zhang ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 148883b2e829SHong Zhang ierr = PetscFree(network->Jv);CHKERRQ(ierr); 14891ad426b7SHong Zhang } 149013c2a604SAdrian Maldonado 149113c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr); 149213c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr); 149313c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr); 149413c2a604SAdrian Maldonado if (network->vertex.sf) { 149513c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr); 149613c2a604SAdrian Maldonado } 149713c2a604SAdrian Maldonado /* edge */ 149813c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr); 149913c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr); 150013c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr); 150113c2a604SAdrian Maldonado if (network->edge.sf) { 150213c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr); 150313c2a604SAdrian Maldonado } 15045f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 15055f2c45f1SShri Abhyankar network->edges = NULL; 15065f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 15075f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 150883b2e829SHong Zhang 15095f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 15105f2c45f1SShri Abhyankar ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 15115f2c45f1SShri Abhyankar ierr = PetscFree(network->header);CHKERRQ(ierr); 15125f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 15135f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15145f2c45f1SShri Abhyankar } 15155f2c45f1SShri Abhyankar 15165f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 15175f2c45f1SShri Abhyankar { 15185f2c45f1SShri Abhyankar PetscErrorCode ierr; 15195f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15205f2c45f1SShri Abhyankar 15215f2c45f1SShri Abhyankar PetscFunctionBegin; 15225f2c45f1SShri Abhyankar ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 15235f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15245f2c45f1SShri Abhyankar } 15255f2c45f1SShri Abhyankar 15265f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 15275f2c45f1SShri Abhyankar { 15285f2c45f1SShri Abhyankar PetscErrorCode ierr; 15295f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15305f2c45f1SShri Abhyankar 15315f2c45f1SShri Abhyankar PetscFunctionBegin; 15325f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 15335f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15345f2c45f1SShri Abhyankar } 15355f2c45f1SShri Abhyankar 15365f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 15375f2c45f1SShri Abhyankar { 15385f2c45f1SShri Abhyankar PetscErrorCode ierr; 15395f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15405f2c45f1SShri Abhyankar 15415f2c45f1SShri Abhyankar PetscFunctionBegin; 15425f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 15435f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15445f2c45f1SShri Abhyankar } 15455f2c45f1SShri Abhyankar 15465f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 15475f2c45f1SShri Abhyankar { 15485f2c45f1SShri Abhyankar PetscErrorCode ierr; 15495f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15505f2c45f1SShri Abhyankar 15515f2c45f1SShri Abhyankar PetscFunctionBegin; 15525f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 15535f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15545f2c45f1SShri Abhyankar } 15555f2c45f1SShri Abhyankar 15565f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 15575f2c45f1SShri Abhyankar { 15585f2c45f1SShri Abhyankar PetscErrorCode ierr; 15595f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15605f2c45f1SShri Abhyankar 15615f2c45f1SShri Abhyankar PetscFunctionBegin; 15625f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 15635f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15645f2c45f1SShri Abhyankar } 1565