1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 25f2c45f1SShri Abhyankar #include <petscdmplex.h> 35f2c45f1SShri Abhyankar #include <petscsf.h> 45f2c45f1SShri Abhyankar 55f2c45f1SShri Abhyankar /*@ 6556ed216SShri Abhyankar DMNetworkGetPlex - Gets the Plex DM associated with this network DM 7556ed216SShri Abhyankar 8556ed216SShri Abhyankar Not collective 9556ed216SShri Abhyankar 10556ed216SShri Abhyankar Input Parameters: 11556ed216SShri Abhyankar + netdm - the dm object 12556ed216SShri Abhyankar - plexmdm - the plex dm object 13556ed216SShri Abhyankar 14556ed216SShri Abhyankar Level: Advanced 15556ed216SShri Abhyankar 16556ed216SShri Abhyankar .seealso: DMNetworkCreate() 17556ed216SShri Abhyankar @*/ 18556ed216SShri Abhyankar PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm) 19556ed216SShri Abhyankar { 20556ed216SShri Abhyankar DM_Network *network = (DM_Network*) netdm->data; 21556ed216SShri Abhyankar 22556ed216SShri Abhyankar PetscFunctionBegin; 23556ed216SShri Abhyankar *plexdm = network->plex; 24556ed216SShri Abhyankar PetscFunctionReturn(0); 25556ed216SShri Abhyankar } 26556ed216SShri Abhyankar 27556ed216SShri Abhyankar /*@ 285f2c45f1SShri Abhyankar DMNetworkSetSizes - Sets the local and global vertices and edges. 295f2c45f1SShri Abhyankar 305f2c45f1SShri Abhyankar Collective on DM 315f2c45f1SShri Abhyankar 325f2c45f1SShri Abhyankar Input Parameters: 335f2c45f1SShri Abhyankar + dm - the dm object 345f2c45f1SShri Abhyankar . nV - number of local vertices 355f2c45f1SShri Abhyankar . nE - number of local edges 365f2c45f1SShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) 375f2c45f1SShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE) 385f2c45f1SShri Abhyankar 395f2c45f1SShri Abhyankar Notes 405f2c45f1SShri Abhyankar If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 415f2c45f1SShri Abhyankar 425f2c45f1SShri Abhyankar You cannot change the sizes once they have been set 435f2c45f1SShri Abhyankar 441b266c99SBarry Smith Level: intermediate 451b266c99SBarry Smith 461b266c99SBarry Smith .seealso: DMNetworkCreate() 475f2c45f1SShri Abhyankar @*/ 485f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt nV, PetscInt nE, PetscInt NV, PetscInt NE) 495f2c45f1SShri Abhyankar { 505f2c45f1SShri Abhyankar PetscErrorCode ierr; 515f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 525f2c45f1SShri Abhyankar PetscInt a[2],b[2]; 535f2c45f1SShri Abhyankar 545f2c45f1SShri Abhyankar PetscFunctionBegin; 555f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 565f2c45f1SShri Abhyankar if (NV > 0) PetscValidLogicalCollectiveInt(dm,NV,4); 575f2c45f1SShri Abhyankar if (NE > 0) PetscValidLogicalCollectiveInt(dm,NE,5); 585f2c45f1SShri Abhyankar if (NV > 0 && nV > NV) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local vertex size %D cannot be larger than global vertex size %D",nV,NV); 595f2c45f1SShri Abhyankar if (NE > 0 && nE > NE) SETERRQ2(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Local edge size %D cannot be larger than global edge size %D",nE,NE); 605f2c45f1SShri Abhyankar if ((network->nNodes >= 0 || network->NNodes >= 0) && (network->nNodes != nV || network->NNodes != NV)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset vertex sizes to %D local %D global after previously setting them to %D local %D global",nV,NV,network->nNodes,network->NNodes); 615f2c45f1SShri Abhyankar if ((network->nEdges >= 0 || network->NEdges >= 0) && (network->nEdges != nE || network->NEdges != NE)) SETERRQ4(PETSC_COMM_SELF,PETSC_ERR_SUP,"Cannot change/reset edge sizes to %D local %D global after previously setting them to %D local %D global",nE,NE,network->nEdges,network->NEdges); 625f2c45f1SShri Abhyankar if (NE < 0 || NV < 0) { 635f2c45f1SShri Abhyankar a[0] = nV; a[1] = nE; 64b2566f29SBarry Smith ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 655f2c45f1SShri Abhyankar NV = b[0]; NE = b[1]; 665f2c45f1SShri Abhyankar } 675f2c45f1SShri Abhyankar network->nNodes = nV; 685f2c45f1SShri Abhyankar network->NNodes = NV; 695f2c45f1SShri Abhyankar network->nEdges = nE; 705f2c45f1SShri Abhyankar network->NEdges = NE; 715f2c45f1SShri Abhyankar PetscFunctionReturn(0); 725f2c45f1SShri Abhyankar } 735f2c45f1SShri Abhyankar 745f2c45f1SShri Abhyankar /*@ 755f2c45f1SShri Abhyankar DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 765f2c45f1SShri Abhyankar 775f2c45f1SShri Abhyankar Logically collective on DM 785f2c45f1SShri Abhyankar 795f2c45f1SShri Abhyankar Input Parameters: 805f2c45f1SShri Abhyankar . edges - list of edges 815f2c45f1SShri Abhyankar 825f2c45f1SShri Abhyankar Notes: 835f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 845f2c45f1SShri Abhyankar not be destroyed before the call to DMNetworkLayoutSetUp 855f2c45f1SShri Abhyankar 865f2c45f1SShri Abhyankar Level: intermediate 875f2c45f1SShri Abhyankar 885f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes 895f2c45f1SShri Abhyankar @*/ 905f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int edgelist[]) 915f2c45f1SShri Abhyankar { 925f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 935f2c45f1SShri Abhyankar 945f2c45f1SShri Abhyankar PetscFunctionBegin; 955f2c45f1SShri Abhyankar network->edges = edgelist; 965f2c45f1SShri Abhyankar PetscFunctionReturn(0); 975f2c45f1SShri Abhyankar } 985f2c45f1SShri Abhyankar 995f2c45f1SShri Abhyankar /*@ 1005f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 1015f2c45f1SShri Abhyankar 1025f2c45f1SShri Abhyankar Collective on DM 1035f2c45f1SShri Abhyankar 1045f2c45f1SShri Abhyankar Input Parameters 1055f2c45f1SShri Abhyankar . DM - the dmnetwork object 1065f2c45f1SShri Abhyankar 1075f2c45f1SShri Abhyankar Notes: 1085f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 1095f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 1105f2c45f1SShri Abhyankar 1115f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 1125f2c45f1SShri Abhyankar 1135f2c45f1SShri Abhyankar Level: intermediate 1145f2c45f1SShri Abhyankar 1155f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 1165f2c45f1SShri Abhyankar @*/ 1175f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 1185f2c45f1SShri Abhyankar { 1195f2c45f1SShri Abhyankar PetscErrorCode ierr; 1205f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1215f2c45f1SShri Abhyankar PetscInt dim = 1; /* One dimensional network */ 1225f2c45f1SShri Abhyankar PetscInt numCorners=2; 1235f2c45f1SShri Abhyankar PetscInt spacedim=2; 1245f2c45f1SShri Abhyankar double *vertexcoords=NULL; 1255f2c45f1SShri Abhyankar PetscInt i; 1265f2c45f1SShri Abhyankar PetscInt ndata; 1275f2c45f1SShri Abhyankar 1285f2c45f1SShri Abhyankar PetscFunctionBegin; 1295f2c45f1SShri Abhyankar if (network->nNodes) { 1309045477aSSatish Balay ierr = PetscCalloc1(numCorners*network->nNodes,&vertexcoords);CHKERRQ(ierr); 1315f2c45f1SShri Abhyankar } 1325f2c45f1SShri Abhyankar ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nNodes,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 1335f2c45f1SShri Abhyankar if (network->nNodes) { 1345f2c45f1SShri Abhyankar ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 1355f2c45f1SShri Abhyankar } 1365f2c45f1SShri Abhyankar ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 1375f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 1385f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 1395f2c45f1SShri Abhyankar 1405f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 1415f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 1425f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 1435f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 1445f2c45f1SShri Abhyankar 14524121865SAdrian Maldonado 14624121865SAdrian Maldonado 1475f2c45f1SShri Abhyankar network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 1486caa05f4SBarry Smith ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 1495f2c45f1SShri Abhyankar for (i = network->pStart; i < network->pEnd; i++) { 1505f2c45f1SShri Abhyankar network->header[i].ndata = 0; 1515f2c45f1SShri Abhyankar ndata = network->header[i].ndata; 1525f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 1535f2c45f1SShri Abhyankar network->header[i].offset[ndata] = 0; 1545f2c45f1SShri Abhyankar } 155854ce69bSBarry Smith ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 1565f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1575f2c45f1SShri Abhyankar } 1585f2c45f1SShri Abhyankar 15994ef8ddeSSatish Balay /*@C 1605f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 1615f2c45f1SShri Abhyankar 1625f2c45f1SShri Abhyankar Logically collective on DM 1635f2c45f1SShri Abhyankar 1645f2c45f1SShri Abhyankar Input Parameters 1655f2c45f1SShri Abhyankar + dm - the network object 1665f2c45f1SShri Abhyankar . name - the component name 1675f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 1685f2c45f1SShri Abhyankar 1695f2c45f1SShri Abhyankar Output Parameters 1705f2c45f1SShri Abhyankar . key - an integer key that defines the component 1715f2c45f1SShri Abhyankar 1725f2c45f1SShri Abhyankar Notes 1735f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 1745f2c45f1SShri Abhyankar 1755f2c45f1SShri Abhyankar Level: intermediate 1765f2c45f1SShri Abhyankar 1775f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 1785f2c45f1SShri Abhyankar @*/ 1795f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 1805f2c45f1SShri Abhyankar { 1815f2c45f1SShri Abhyankar PetscErrorCode ierr; 1825f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1835f2c45f1SShri Abhyankar DMNetworkComponent *component=&network->component[network->ncomponent]; 1845f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 1855f2c45f1SShri Abhyankar PetscInt i; 1865f2c45f1SShri Abhyankar 1875f2c45f1SShri Abhyankar PetscFunctionBegin; 1885f2c45f1SShri Abhyankar 1895f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 1905f2c45f1SShri Abhyankar ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 1915f2c45f1SShri Abhyankar if (flg) { 1925f2c45f1SShri Abhyankar *key = i; 1935f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1945f2c45f1SShri Abhyankar } 1955f2c45f1SShri Abhyankar } 1965f2c45f1SShri Abhyankar 1975f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 1985f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 1995f2c45f1SShri Abhyankar *key = network->ncomponent; 2005f2c45f1SShri Abhyankar network->ncomponent++; 2015f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2025f2c45f1SShri Abhyankar } 2035f2c45f1SShri Abhyankar 2045f2c45f1SShri Abhyankar /*@ 2055f2c45f1SShri Abhyankar DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 2065f2c45f1SShri Abhyankar 2075f2c45f1SShri Abhyankar Not Collective 2085f2c45f1SShri Abhyankar 2095f2c45f1SShri Abhyankar Input Parameters: 2105f2c45f1SShri Abhyankar + dm - The DMNetwork object 2115f2c45f1SShri Abhyankar 2125f2c45f1SShri Abhyankar Output Paramters: 2135f2c45f1SShri Abhyankar + vStart - The first vertex point 2145f2c45f1SShri Abhyankar - vEnd - One beyond the last vertex point 2155f2c45f1SShri Abhyankar 2165f2c45f1SShri Abhyankar Level: intermediate 2175f2c45f1SShri Abhyankar 2185f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange 2195f2c45f1SShri Abhyankar @*/ 2205f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 2215f2c45f1SShri Abhyankar { 2225f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2235f2c45f1SShri Abhyankar 2245f2c45f1SShri Abhyankar PetscFunctionBegin; 2255f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 2265f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 2275f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2285f2c45f1SShri Abhyankar } 2295f2c45f1SShri Abhyankar 2305f2c45f1SShri Abhyankar /*@ 2315f2c45f1SShri Abhyankar DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 2325f2c45f1SShri Abhyankar 2335f2c45f1SShri Abhyankar Not Collective 2345f2c45f1SShri Abhyankar 2355f2c45f1SShri Abhyankar Input Parameters: 2365f2c45f1SShri Abhyankar + dm - The DMNetwork object 2375f2c45f1SShri Abhyankar 2385f2c45f1SShri Abhyankar Output Paramters: 2395f2c45f1SShri Abhyankar + eStart - The first edge point 2405f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 2415f2c45f1SShri Abhyankar 2425f2c45f1SShri Abhyankar Level: intermediate 2435f2c45f1SShri Abhyankar 2445f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange 2455f2c45f1SShri Abhyankar @*/ 2465f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 2475f2c45f1SShri Abhyankar { 2485f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2495f2c45f1SShri Abhyankar 2505f2c45f1SShri Abhyankar PetscFunctionBegin; 2515f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 2525f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 2535f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2545f2c45f1SShri Abhyankar } 2555f2c45f1SShri Abhyankar 2565f2c45f1SShri Abhyankar /*@ 2575f2c45f1SShri Abhyankar DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 2585f2c45f1SShri Abhyankar 2595f2c45f1SShri Abhyankar Not Collective 2605f2c45f1SShri Abhyankar 2615f2c45f1SShri Abhyankar Input Parameters: 2625f2c45f1SShri Abhyankar + dm - The DMNetwork object 2635f2c45f1SShri Abhyankar . p - vertex/edge point 2645f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component 2655f2c45f1SShri Abhyankar - compvalue - pointer to the data structure for the component 2665f2c45f1SShri Abhyankar 2675f2c45f1SShri Abhyankar Level: intermediate 2685f2c45f1SShri Abhyankar 2695f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 2705f2c45f1SShri Abhyankar @*/ 2715f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 2725f2c45f1SShri Abhyankar { 2735f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 27443a39a44SBarry Smith DMNetworkComponent *component = &network->component[componentkey]; 2755f2c45f1SShri Abhyankar DMNetworkComponentHeader header = &network->header[p]; 2765f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue = &network->cvalue[p]; 2775f2c45f1SShri Abhyankar PetscErrorCode ierr; 2785f2c45f1SShri Abhyankar 2795f2c45f1SShri Abhyankar PetscFunctionBegin; 280fa58f0a9SHong 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); 281fa58f0a9SHong Zhang 28243a39a44SBarry Smith header->size[header->ndata] = component->size; 28343a39a44SBarry Smith ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 2845f2c45f1SShri Abhyankar header->key[header->ndata] = componentkey; 2855f2c45f1SShri Abhyankar if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 2865f2c45f1SShri Abhyankar 2875f2c45f1SShri Abhyankar cvalue->data[header->ndata] = (void*)compvalue; 2885f2c45f1SShri Abhyankar header->ndata++; 2895f2c45f1SShri Abhyankar PetscFunctionReturn(0); 2905f2c45f1SShri Abhyankar } 2915f2c45f1SShri Abhyankar 2925f2c45f1SShri Abhyankar /*@ 2935f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 2945f2c45f1SShri Abhyankar 2955f2c45f1SShri Abhyankar Not Collective 2965f2c45f1SShri Abhyankar 2975f2c45f1SShri Abhyankar Input Parameters: 2985f2c45f1SShri Abhyankar + dm - The DMNetwork object 2995f2c45f1SShri Abhyankar . p - vertex/edge point 3005f2c45f1SShri Abhyankar 3015f2c45f1SShri Abhyankar Output Parameters: 3025f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 3035f2c45f1SShri Abhyankar 3045f2c45f1SShri Abhyankar Level: intermediate 3055f2c45f1SShri Abhyankar 3065f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 3075f2c45f1SShri Abhyankar @*/ 3085f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 3095f2c45f1SShri Abhyankar { 3105f2c45f1SShri Abhyankar PetscErrorCode ierr; 3115f2c45f1SShri Abhyankar PetscInt offset; 3125f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3135f2c45f1SShri Abhyankar 3145f2c45f1SShri Abhyankar PetscFunctionBegin; 3155f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 3165f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 3175f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3185f2c45f1SShri Abhyankar } 3195f2c45f1SShri Abhyankar 3205f2c45f1SShri Abhyankar /*@ 3215f2c45f1SShri Abhyankar DMNetworkGetComponentTypeOffset - Gets the type along with the offset for indexing the 3225f2c45f1SShri Abhyankar component value from the component data array 3235f2c45f1SShri Abhyankar 3245f2c45f1SShri Abhyankar Not Collective 3255f2c45f1SShri Abhyankar 3265f2c45f1SShri Abhyankar Input Parameters: 3275f2c45f1SShri Abhyankar + dm - The DMNetwork object 3285f2c45f1SShri Abhyankar . p - vertex/edge point 3295f2c45f1SShri Abhyankar - compnum - component number 3305f2c45f1SShri Abhyankar 3315f2c45f1SShri Abhyankar Output Parameters: 3325f2c45f1SShri Abhyankar + compkey - the key obtained when registering the component 3335f2c45f1SShri Abhyankar - offset - offset into the component data array associated with the vertex/edge point 3345f2c45f1SShri Abhyankar 3355f2c45f1SShri Abhyankar Notes: 3365f2c45f1SShri Abhyankar Typical usage: 3375f2c45f1SShri Abhyankar 3385f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray(dm, &arr); 3395f2c45f1SShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 3405f2c45f1SShri Abhyankar Loop over vertices or edges 3415f2c45f1SShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 3425f2c45f1SShri Abhyankar Loop over numcomps 3435f2c45f1SShri Abhyankar DMNetworkGetComponentTypeOffset(dm,v,compnum,&key,&offset); 3445f2c45f1SShri Abhyankar compdata = (UserCompDataType)(arr+offset); 3455f2c45f1SShri Abhyankar 3465f2c45f1SShri Abhyankar Level: intermediate 3475f2c45f1SShri Abhyankar 3485f2c45f1SShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 3495f2c45f1SShri Abhyankar @*/ 3505f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentTypeOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 3515f2c45f1SShri Abhyankar { 3525f2c45f1SShri Abhyankar PetscErrorCode ierr; 3535f2c45f1SShri Abhyankar PetscInt offsetp; 3545f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 3555f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3565f2c45f1SShri Abhyankar 3575f2c45f1SShri Abhyankar PetscFunctionBegin; 3585f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 3595f2c45f1SShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 360d36f4e81SHong Zhang if (compkey) *compkey = header->key[compnum]; 361d36f4e81SHong Zhang if (offset) *offset = offsetp+network->dataheadersize+header->offset[compnum]; 3625f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3635f2c45f1SShri Abhyankar } 3645f2c45f1SShri Abhyankar 3655f2c45f1SShri Abhyankar /*@ 3665f2c45f1SShri Abhyankar DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 3675f2c45f1SShri Abhyankar 3685f2c45f1SShri Abhyankar Not Collective 3695f2c45f1SShri Abhyankar 3705f2c45f1SShri Abhyankar Input Parameters: 3715f2c45f1SShri Abhyankar + dm - The DMNetwork object 3725f2c45f1SShri Abhyankar - p - the edge/vertex point 3735f2c45f1SShri Abhyankar 3745f2c45f1SShri Abhyankar Output Parameters: 3755f2c45f1SShri Abhyankar . offset - the offset 3765f2c45f1SShri Abhyankar 3775f2c45f1SShri Abhyankar Level: intermediate 3785f2c45f1SShri Abhyankar 3795f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 3805f2c45f1SShri Abhyankar @*/ 3815f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 3825f2c45f1SShri Abhyankar { 3835f2c45f1SShri Abhyankar PetscErrorCode ierr; 3845f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 3855f2c45f1SShri Abhyankar 3865f2c45f1SShri Abhyankar PetscFunctionBegin; 3875f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr); 3885f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3895f2c45f1SShri Abhyankar } 3905f2c45f1SShri Abhyankar 3915f2c45f1SShri Abhyankar /*@ 3925f2c45f1SShri Abhyankar DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 3935f2c45f1SShri Abhyankar 3945f2c45f1SShri Abhyankar Not Collective 3955f2c45f1SShri Abhyankar 3965f2c45f1SShri Abhyankar Input Parameters: 3975f2c45f1SShri Abhyankar + dm - The DMNetwork object 3985f2c45f1SShri Abhyankar - p - the edge/vertex point 3995f2c45f1SShri Abhyankar 4005f2c45f1SShri Abhyankar Output Parameters: 4015f2c45f1SShri Abhyankar . offsetg - the offset 4025f2c45f1SShri Abhyankar 4035f2c45f1SShri Abhyankar Level: intermediate 4045f2c45f1SShri Abhyankar 4055f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 4065f2c45f1SShri Abhyankar @*/ 4075f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 4085f2c45f1SShri Abhyankar { 4095f2c45f1SShri Abhyankar PetscErrorCode ierr; 4105f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4115f2c45f1SShri Abhyankar 4125f2c45f1SShri Abhyankar PetscFunctionBegin; 4135f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr); 414dd6f46cdSHong Zhang if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost node */ 4155f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4165f2c45f1SShri Abhyankar } 4175f2c45f1SShri Abhyankar 41824121865SAdrian Maldonado /*@ 41924121865SAdrian Maldonado DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector. 42024121865SAdrian Maldonado 42124121865SAdrian Maldonado Not Collective 42224121865SAdrian Maldonado 42324121865SAdrian Maldonado Input Parameters: 42424121865SAdrian Maldonado + dm - The DMNetwork object 42524121865SAdrian Maldonado - p - the edge point 42624121865SAdrian Maldonado 42724121865SAdrian Maldonado Output Parameters: 42824121865SAdrian Maldonado . offset - the offset 42924121865SAdrian Maldonado 43024121865SAdrian Maldonado Level: intermediate 43124121865SAdrian Maldonado 43224121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 43324121865SAdrian Maldonado @*/ 43424121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset) 43524121865SAdrian Maldonado { 43624121865SAdrian Maldonado PetscErrorCode ierr; 43724121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 43824121865SAdrian Maldonado 43924121865SAdrian Maldonado PetscFunctionBegin; 44024121865SAdrian Maldonado 44124121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr); 44224121865SAdrian Maldonado PetscFunctionReturn(0); 44324121865SAdrian Maldonado } 44424121865SAdrian Maldonado 44524121865SAdrian Maldonado /*@ 44624121865SAdrian Maldonado DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector. 44724121865SAdrian Maldonado 44824121865SAdrian Maldonado Not Collective 44924121865SAdrian Maldonado 45024121865SAdrian Maldonado Input Parameters: 45124121865SAdrian Maldonado + dm - The DMNetwork object 45224121865SAdrian Maldonado - p - the vertex point 45324121865SAdrian Maldonado 45424121865SAdrian Maldonado Output Parameters: 45524121865SAdrian Maldonado . offset - the offset 45624121865SAdrian Maldonado 45724121865SAdrian Maldonado Level: intermediate 45824121865SAdrian Maldonado 45924121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 46024121865SAdrian Maldonado @*/ 46124121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset) 46224121865SAdrian Maldonado { 46324121865SAdrian Maldonado PetscErrorCode ierr; 46424121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 46524121865SAdrian Maldonado 46624121865SAdrian Maldonado PetscFunctionBegin; 46724121865SAdrian Maldonado 46824121865SAdrian Maldonado p -= network->vStart; 46924121865SAdrian Maldonado 47024121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr); 47124121865SAdrian Maldonado PetscFunctionReturn(0); 47224121865SAdrian Maldonado } 4735f2c45f1SShri Abhyankar /*@ 4745f2c45f1SShri Abhyankar DMNetworkAddNumVariables - Add number of variables associated with a given point. 4755f2c45f1SShri Abhyankar 4765f2c45f1SShri Abhyankar Not Collective 4775f2c45f1SShri Abhyankar 4785f2c45f1SShri Abhyankar Input Parameters: 4795f2c45f1SShri Abhyankar + dm - The DMNetworkObject 4805f2c45f1SShri Abhyankar . p - the vertex/edge point 4815f2c45f1SShri Abhyankar - nvar - number of additional variables 4825f2c45f1SShri Abhyankar 4835f2c45f1SShri Abhyankar Level: intermediate 4845f2c45f1SShri Abhyankar 4855f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables 4865f2c45f1SShri Abhyankar @*/ 4875f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 4885f2c45f1SShri Abhyankar { 4895f2c45f1SShri Abhyankar PetscErrorCode ierr; 4905f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4915f2c45f1SShri Abhyankar 4925f2c45f1SShri Abhyankar PetscFunctionBegin; 4935f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 4945f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4955f2c45f1SShri Abhyankar } 4965f2c45f1SShri Abhyankar 49727f51fceSHong Zhang /*@ 49827f51fceSHong Zhang DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point. 49927f51fceSHong Zhang 50027f51fceSHong Zhang Not Collective 50127f51fceSHong Zhang 50227f51fceSHong Zhang Input Parameters: 50327f51fceSHong Zhang + dm - The DMNetworkObject 50427f51fceSHong Zhang - p - the vertex/edge point 50527f51fceSHong Zhang 50627f51fceSHong Zhang Output Parameters: 50727f51fceSHong Zhang . nvar - number of variables 50827f51fceSHong Zhang 50927f51fceSHong Zhang Level: intermediate 51027f51fceSHong Zhang 51127f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables 51227f51fceSHong Zhang @*/ 51327f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar) 51427f51fceSHong Zhang { 51527f51fceSHong Zhang PetscErrorCode ierr; 51627f51fceSHong Zhang DM_Network *network = (DM_Network*)dm->data; 51727f51fceSHong Zhang 51827f51fceSHong Zhang PetscFunctionBegin; 51927f51fceSHong Zhang ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 52027f51fceSHong Zhang PetscFunctionReturn(0); 52127f51fceSHong Zhang } 52227f51fceSHong Zhang 5235f2c45f1SShri Abhyankar /*@ 5245f2c45f1SShri Abhyankar DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 5255f2c45f1SShri Abhyankar 5265f2c45f1SShri Abhyankar Not Collective 5275f2c45f1SShri Abhyankar 5285f2c45f1SShri Abhyankar Input Parameters: 5295f2c45f1SShri Abhyankar + dm - The DMNetworkObject 5305f2c45f1SShri Abhyankar . p - the vertex/edge point 5315f2c45f1SShri Abhyankar - nvar - number of variables 5325f2c45f1SShri Abhyankar 5335f2c45f1SShri Abhyankar Level: intermediate 5345f2c45f1SShri Abhyankar 5355f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables 5365f2c45f1SShri Abhyankar @*/ 5375f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 5385f2c45f1SShri Abhyankar { 5395f2c45f1SShri Abhyankar PetscErrorCode ierr; 5405f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5415f2c45f1SShri Abhyankar 5425f2c45f1SShri Abhyankar PetscFunctionBegin; 5435f2c45f1SShri Abhyankar ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 5445f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5455f2c45f1SShri Abhyankar } 5465f2c45f1SShri Abhyankar 5475f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This 5485f2c45f1SShri Abhyankar function is called during DMSetUp() */ 5495f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 5505f2c45f1SShri Abhyankar { 5515f2c45f1SShri Abhyankar PetscErrorCode ierr; 5525f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5535f2c45f1SShri Abhyankar PetscInt arr_size; 5545f2c45f1SShri Abhyankar PetscInt p,offset,offsetp; 5555f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 5565f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 5575f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 5585f2c45f1SShri Abhyankar PetscInt ncomp, i; 5595f2c45f1SShri Abhyankar 5605f2c45f1SShri Abhyankar PetscFunctionBegin; 5615f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 5625f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 56375b160a0SShri Abhyankar ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 5645f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 5655f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 5665f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 5675f2c45f1SShri Abhyankar /* Copy header */ 5685f2c45f1SShri Abhyankar header = &network->header[p]; 569302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 5705f2c45f1SShri Abhyankar /* Copy data */ 5715f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 5725f2c45f1SShri Abhyankar ncomp = header->ndata; 5735f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 5745f2c45f1SShri Abhyankar offset = offsetp + network->dataheadersize + header->offset[i]; 575302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 5765f2c45f1SShri Abhyankar } 5775f2c45f1SShri Abhyankar } 5785f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5795f2c45f1SShri Abhyankar } 5805f2c45f1SShri Abhyankar 5815f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 5825f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm) 5835f2c45f1SShri Abhyankar { 5845f2c45f1SShri Abhyankar PetscErrorCode ierr; 5855f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5865f2c45f1SShri Abhyankar 5875f2c45f1SShri Abhyankar PetscFunctionBegin; 5885f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 5895f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5905f2c45f1SShri Abhyankar } 5915f2c45f1SShri Abhyankar 5925f2c45f1SShri Abhyankar /*@C 5935f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray - Returns the component data array 5945f2c45f1SShri Abhyankar 5955f2c45f1SShri Abhyankar Not Collective 5965f2c45f1SShri Abhyankar 5975f2c45f1SShri Abhyankar Input Parameters: 5985f2c45f1SShri Abhyankar . dm - The DMNetwork Object 5995f2c45f1SShri Abhyankar 6005f2c45f1SShri Abhyankar Output Parameters: 6015f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components 6025f2c45f1SShri Abhyankar 6035f2c45f1SShri Abhyankar Level: intermediate 6045f2c45f1SShri Abhyankar 6055f2c45f1SShri Abhyankar .seealso: DMNetworkGetComponentTypeOffset, DMNetworkGetNumComponents 6065f2c45f1SShri Abhyankar @*/ 6075f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 6085f2c45f1SShri Abhyankar { 6095f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 6105f2c45f1SShri Abhyankar 6115f2c45f1SShri Abhyankar PetscFunctionBegin; 6125f2c45f1SShri Abhyankar *componentdataarray = network->componentdataarray; 6135f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6145f2c45f1SShri Abhyankar } 6155f2c45f1SShri Abhyankar 61624121865SAdrian Maldonado /* Get a subsection from a range of points */ 61724121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection) 61824121865SAdrian Maldonado { 61924121865SAdrian Maldonado PetscErrorCode ierr; 62024121865SAdrian Maldonado PetscInt i, nvar; 62124121865SAdrian Maldonado 62224121865SAdrian Maldonado PetscFunctionBegin; 62324121865SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr); 62424121865SAdrian Maldonado ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr); 62524121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 62624121865SAdrian Maldonado ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr); 62724121865SAdrian Maldonado ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr); 62824121865SAdrian Maldonado } 62924121865SAdrian Maldonado 63024121865SAdrian Maldonado ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr); 63124121865SAdrian Maldonado PetscFunctionReturn(0); 63224121865SAdrian Maldonado } 63324121865SAdrian Maldonado 63424121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */ 63524121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 63624121865SAdrian Maldonado { 63724121865SAdrian Maldonado PetscErrorCode ierr; 63824121865SAdrian Maldonado PetscInt i, *subpoints; 63924121865SAdrian Maldonado 64024121865SAdrian Maldonado PetscFunctionBegin; 64124121865SAdrian Maldonado /* Create index sets to map from "points" to "subpoints" */ 64224121865SAdrian Maldonado ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr); 64324121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 64424121865SAdrian Maldonado subpoints[i - pstart] = i; 64524121865SAdrian Maldonado } 646459726d8SSatish Balay ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr); 64724121865SAdrian Maldonado ierr = PetscFree(subpoints);CHKERRQ(ierr); 64824121865SAdrian Maldonado PetscFunctionReturn(0); 64924121865SAdrian Maldonado } 65024121865SAdrian Maldonado 65124121865SAdrian Maldonado /*@ 65224121865SAdrian Maldonado DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute. 65324121865SAdrian Maldonado 65424121865SAdrian Maldonado Collective 65524121865SAdrian Maldonado 65624121865SAdrian Maldonado Input Parameters: 65724121865SAdrian Maldonado . dm - The DMNetworkObject 65824121865SAdrian Maldonado 65924121865SAdrian Maldonado Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 66024121865SAdrian Maldonado 66124121865SAdrian Maldonado points = [0 1 2 3 4 5 6] 66224121865SAdrian Maldonado 66324121865SAdrian Maldonado where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]). 66424121865SAdrian Maldonado 66524121865SAdrian Maldonado With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 66624121865SAdrian Maldonado 66724121865SAdrian Maldonado Level: intermediate 66824121865SAdrian Maldonado 66924121865SAdrian Maldonado @*/ 67024121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 67124121865SAdrian Maldonado { 67224121865SAdrian Maldonado PetscErrorCode ierr; 67324121865SAdrian Maldonado MPI_Comm comm; 6749852e123SBarry Smith PetscMPIInt rank, size; 67524121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 67624121865SAdrian Maldonado 677eab1376dSHong Zhang PetscFunctionBegin; 67824121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 67924121865SAdrian Maldonado ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 6809852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 68124121865SAdrian Maldonado 68224121865SAdrian Maldonado /* Create maps for vertices and edges */ 68324121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 68424121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr); 68524121865SAdrian Maldonado 68624121865SAdrian Maldonado /* Create local sub-sections */ 68724121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr); 68824121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr); 68924121865SAdrian Maldonado 6909852e123SBarry Smith if (size > 1) { 69124121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 69224121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr); 69324121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr); 69424121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr); 69524121865SAdrian Maldonado } else { 69624121865SAdrian Maldonado /* create structures for vertex */ 69724121865SAdrian Maldonado ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr); 69824121865SAdrian Maldonado /* create structures for edge */ 69924121865SAdrian Maldonado ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr); 70024121865SAdrian Maldonado } 70124121865SAdrian Maldonado 70224121865SAdrian Maldonado 70324121865SAdrian Maldonado /* Add viewers */ 70424121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr); 70524121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr); 70624121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr); 70724121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr); 70824121865SAdrian Maldonado 70924121865SAdrian Maldonado PetscFunctionReturn(0); 71024121865SAdrian Maldonado } 7115f2c45f1SShri Abhyankar /*@ 7125f2c45f1SShri Abhyankar DMNetworkDistribute - Distributes the network and moves associated component data. 7135f2c45f1SShri Abhyankar 7145f2c45f1SShri Abhyankar Collective 7155f2c45f1SShri Abhyankar 7165f2c45f1SShri Abhyankar Input Parameter: 717*d3464fd4SAdrian Maldonado + DM - the DMNetwork object 7185f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default 7195f2c45f1SShri Abhyankar 7205f2c45f1SShri Abhyankar Notes: 7218b171c8eSHong Zhang Distributes the network with <overlap>-overlapping partitioning of the edges. 7225f2c45f1SShri Abhyankar 7235f2c45f1SShri Abhyankar Level: intermediate 7245f2c45f1SShri Abhyankar 7255f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 7265f2c45f1SShri Abhyankar @*/ 727*d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm, PetscInt overlap) 7285f2c45f1SShri Abhyankar { 729*d3464fd4SAdrian Maldonado MPI_Comm comm; 7305f2c45f1SShri Abhyankar PetscErrorCode ierr; 731*d3464fd4SAdrian Maldonado PetscMPIInt size; 732*d3464fd4SAdrian Maldonado DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data); 733*d3464fd4SAdrian Maldonado DM_Network *newDMnetwork; 7345f2c45f1SShri Abhyankar PetscSF pointsf; 7355f2c45f1SShri Abhyankar DM newDM; 7365f2c45f1SShri Abhyankar 7375f2c45f1SShri Abhyankar PetscFunctionBegin; 738*d3464fd4SAdrian Maldonado 739*d3464fd4SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr); 740*d3464fd4SAdrian Maldonado ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 741*d3464fd4SAdrian Maldonado 742*d3464fd4SAdrian Maldonado if (size == 1) PetscFunctionReturn(0); 743*d3464fd4SAdrian Maldonado 744*d3464fd4SAdrian Maldonado ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr); 7455f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 7465f2c45f1SShri Abhyankar newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 7475f2c45f1SShri Abhyankar /* Distribute plex dm and dof section */ 74880cf41d5SMatthew G. Knepley ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 7495f2c45f1SShri Abhyankar /* Distribute dof section */ 750*d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr); 7515f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 752*d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr); 7535f2c45f1SShri Abhyankar /* Distribute data and associated section */ 75431da1fc8SHong Zhang ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 75524121865SAdrian Maldonado 7565f2c45f1SShri Abhyankar 7575f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 7585f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 7595f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 7605f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 7615f2c45f1SShri Abhyankar newDMnetwork->nNodes = newDMnetwork->vEnd - newDMnetwork->vStart; 7625f2c45f1SShri Abhyankar newDMnetwork->NNodes = oldDMnetwork->NNodes; 7635f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 76424121865SAdrian Maldonado 7655f2c45f1SShri Abhyankar /* Set Dof section as the default section for dm */ 7665f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 7675f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 7685f2c45f1SShri Abhyankar 76924121865SAdrian Maldonado /* Destroy point SF */ 77024121865SAdrian Maldonado ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 77124121865SAdrian Maldonado 772*d3464fd4SAdrian Maldonado ierr = DMDestroy(dm);CHKERRQ(ierr); 773*d3464fd4SAdrian Maldonado *dm = newDM; 774*d3464fd4SAdrian Maldonado 7755f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7765f2c45f1SShri Abhyankar } 7775f2c45f1SShri Abhyankar 77824121865SAdrian Maldonado /*@C 77924121865SAdrian Maldonado PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering. 78024121865SAdrian Maldonado 78124121865SAdrian Maldonado Input Parameters: 78224121865SAdrian Maldonado + masterSF - the original SF structure 78324121865SAdrian Maldonado - map - a ISLocalToGlobal mapping that contains the subset of points 78424121865SAdrian Maldonado 78524121865SAdrian Maldonado Output Parameters: 78624121865SAdrian Maldonado . subSF - a subset of the masterSF for the desired subset. 78724121865SAdrian Maldonado */ 78824121865SAdrian Maldonado 78924121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) { 79024121865SAdrian Maldonado 79124121865SAdrian Maldonado PetscErrorCode ierr; 79224121865SAdrian Maldonado PetscInt nroots, nleaves, *ilocal_sub; 79324121865SAdrian Maldonado PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 79424121865SAdrian Maldonado PetscInt *local_points, *remote_points; 79524121865SAdrian Maldonado PetscSFNode *iremote_sub; 79624121865SAdrian Maldonado const PetscInt *ilocal; 79724121865SAdrian Maldonado const PetscSFNode *iremote; 79824121865SAdrian Maldonado 79924121865SAdrian Maldonado PetscFunctionBegin; 80024121865SAdrian Maldonado ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 80124121865SAdrian Maldonado 80224121865SAdrian Maldonado /* Look for leaves that pertain to the subset of points. Get the local ordering */ 80324121865SAdrian Maldonado ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr); 80424121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr); 80524121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 80624121865SAdrian Maldonado if (ilocal_map[i] != -1) nleaves_sub += 1; 80724121865SAdrian Maldonado } 80824121865SAdrian Maldonado /* Re-number ilocal with subset numbering. Need information from roots */ 80924121865SAdrian Maldonado ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr); 81024121865SAdrian Maldonado for (i = 0; i < nroots; i++) local_points[i] = i; 81124121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr); 81224121865SAdrian Maldonado ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 81324121865SAdrian Maldonado ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 81424121865SAdrian Maldonado /* Fill up graph using local (that is, local to the subset) numbering. */ 8154b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr); 8164b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr); 81724121865SAdrian Maldonado nleaves_sub = 0; 81824121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 81924121865SAdrian Maldonado if (ilocal_map[i] != -1) { 82024121865SAdrian Maldonado ilocal_sub[nleaves_sub] = ilocal_map[i]; 8214b70a8deSAdrian Maldonado iremote_sub[nleaves_sub].rank = iremote[i].rank; 82224121865SAdrian Maldonado iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 82324121865SAdrian Maldonado nleaves_sub += 1; 82424121865SAdrian Maldonado } 82524121865SAdrian Maldonado } 82624121865SAdrian Maldonado ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr); 82724121865SAdrian Maldonado ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr); 82824121865SAdrian Maldonado 82924121865SAdrian Maldonado /* Create new subSF */ 83024121865SAdrian Maldonado ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr); 83124121865SAdrian Maldonado ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr); 8324b70a8deSAdrian Maldonado ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr); 83324121865SAdrian Maldonado ierr = PetscFree(ilocal_map);CHKERRQ(ierr); 8344b70a8deSAdrian Maldonado ierr = PetscFree(iremote_sub);CHKERRQ(ierr); 83524121865SAdrian Maldonado 83624121865SAdrian Maldonado PetscFunctionReturn(0); 83724121865SAdrian Maldonado } 83824121865SAdrian Maldonado 83924121865SAdrian Maldonado 8405f2c45f1SShri Abhyankar /*@C 8415f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 8425f2c45f1SShri Abhyankar 8435f2c45f1SShri Abhyankar Not Collective 8445f2c45f1SShri Abhyankar 8455f2c45f1SShri Abhyankar Input Parameters: 8465f2c45f1SShri Abhyankar + dm - The DMNetwork object 8475f2c45f1SShri Abhyankar - p - the vertex point 8485f2c45f1SShri Abhyankar 8495f2c45f1SShri Abhyankar Output Paramters: 8505f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 8515f2c45f1SShri Abhyankar - edges - List of edge points 8525f2c45f1SShri Abhyankar 8535f2c45f1SShri Abhyankar Level: intermediate 8545f2c45f1SShri Abhyankar 8555f2c45f1SShri Abhyankar Fortran Notes: 8565f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 8575f2c45f1SShri Abhyankar include petsc.h90 in your code. 8585f2c45f1SShri Abhyankar 8595f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes 8605f2c45f1SShri Abhyankar @*/ 8615f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 8625f2c45f1SShri Abhyankar { 8635f2c45f1SShri Abhyankar PetscErrorCode ierr; 8645f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8655f2c45f1SShri Abhyankar 8665f2c45f1SShri Abhyankar PetscFunctionBegin; 8675f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 8685f2c45f1SShri Abhyankar ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 8695f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8705f2c45f1SShri Abhyankar } 8715f2c45f1SShri Abhyankar 8725f2c45f1SShri Abhyankar /*@C 873596e729fSHong Zhang DMNetworkGetConnectedNodes - Return the connected vertices for this edge point 8745f2c45f1SShri Abhyankar 8755f2c45f1SShri Abhyankar Not Collective 8765f2c45f1SShri Abhyankar 8775f2c45f1SShri Abhyankar Input Parameters: 8785f2c45f1SShri Abhyankar + dm - The DMNetwork object 8795f2c45f1SShri Abhyankar - p - the edge point 8805f2c45f1SShri Abhyankar 8815f2c45f1SShri Abhyankar Output Paramters: 8825f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 8835f2c45f1SShri Abhyankar 8845f2c45f1SShri Abhyankar Level: intermediate 8855f2c45f1SShri Abhyankar 8865f2c45f1SShri Abhyankar Fortran Notes: 8875f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 8885f2c45f1SShri Abhyankar include petsc.h90 in your code. 8895f2c45f1SShri Abhyankar 8905f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 8915f2c45f1SShri Abhyankar @*/ 8925f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetConnectedNodes(DM dm,PetscInt edge,const PetscInt *vertices[]) 8935f2c45f1SShri Abhyankar { 8945f2c45f1SShri Abhyankar PetscErrorCode ierr; 8955f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8965f2c45f1SShri Abhyankar 8975f2c45f1SShri Abhyankar PetscFunctionBegin; 8985f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 8995f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9005f2c45f1SShri Abhyankar } 9015f2c45f1SShri Abhyankar 9025f2c45f1SShri Abhyankar /*@ 9035f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 9045f2c45f1SShri Abhyankar 9055f2c45f1SShri Abhyankar Not Collective 9065f2c45f1SShri Abhyankar 9075f2c45f1SShri Abhyankar Input Parameters: 9085f2c45f1SShri Abhyankar + dm - The DMNetwork object 9095f2c45f1SShri Abhyankar . p - the vertex point 9105f2c45f1SShri Abhyankar 9115f2c45f1SShri Abhyankar Output Parameter: 9125f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 9135f2c45f1SShri Abhyankar 9145f2c45f1SShri Abhyankar Level: intermediate 9155f2c45f1SShri Abhyankar 9165f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetConnectedNodes, DMNetworkGetVertexRange 9175f2c45f1SShri Abhyankar @*/ 9185f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 9195f2c45f1SShri Abhyankar { 9205f2c45f1SShri Abhyankar PetscErrorCode ierr; 9215f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9225f2c45f1SShri Abhyankar PetscInt offsetg; 9235f2c45f1SShri Abhyankar PetscSection sectiong; 9245f2c45f1SShri Abhyankar 9255f2c45f1SShri Abhyankar PetscFunctionBegin; 9265f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 9275f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 9285f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 9295f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 9305f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9315f2c45f1SShri Abhyankar } 9325f2c45f1SShri Abhyankar 9335f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 9345f2c45f1SShri Abhyankar { 9355f2c45f1SShri Abhyankar PetscErrorCode ierr; 9365f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 9375f2c45f1SShri Abhyankar 9385f2c45f1SShri Abhyankar PetscFunctionBegin; 9395f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 9405f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 9415f2c45f1SShri Abhyankar 9425f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 9435f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 9445f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9455f2c45f1SShri Abhyankar } 9465f2c45f1SShri Abhyankar 9471ad426b7SHong Zhang /*@ 94817df6e9eSHong Zhang DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 9491ad426b7SHong Zhang -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 9501ad426b7SHong Zhang 9511ad426b7SHong Zhang Collective 9521ad426b7SHong Zhang 9531ad426b7SHong Zhang Input Parameters: 95483b2e829SHong Zhang + dm - The DMNetwork object 95583b2e829SHong Zhang . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 95683b2e829SHong Zhang - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 9571ad426b7SHong Zhang 9581ad426b7SHong Zhang Level: intermediate 9591ad426b7SHong Zhang 9601ad426b7SHong Zhang @*/ 96183b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 9621ad426b7SHong Zhang { 9631ad426b7SHong Zhang DM_Network *network=(DM_Network*)dm->data; 9641ad426b7SHong Zhang 9651ad426b7SHong Zhang PetscFunctionBegin; 96683b2e829SHong Zhang network->userEdgeJacobian = eflg; 96783b2e829SHong Zhang network->userVertexJacobian = vflg; 9681ad426b7SHong Zhang PetscFunctionReturn(0); 9691ad426b7SHong Zhang } 9701ad426b7SHong Zhang 9711ad426b7SHong Zhang /*@ 97283b2e829SHong Zhang DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 97383b2e829SHong Zhang 97483b2e829SHong Zhang Not Collective 97583b2e829SHong Zhang 97683b2e829SHong Zhang Input Parameters: 97783b2e829SHong Zhang + dm - The DMNetwork object 97883b2e829SHong Zhang . p - the edge point 9793e97b6e8SHong Zhang - J - array (size = 3) of Jacobian submatrices for this edge point: 9803e97b6e8SHong Zhang J[0]: this edge 9813e97b6e8SHong Zhang J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedNodes() 98283b2e829SHong Zhang 98383b2e829SHong Zhang Level: intermediate 98483b2e829SHong Zhang 98583b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix 98683b2e829SHong Zhang @*/ 98783b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 98883b2e829SHong Zhang { 98983b2e829SHong Zhang PetscErrorCode ierr; 99083b2e829SHong Zhang DM_Network *network=(DM_Network*)dm->data; 99183b2e829SHong Zhang 99283b2e829SHong Zhang PetscFunctionBegin; 993883e35e8SHong Zhang if (!network->userEdgeJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 99483b2e829SHong Zhang if (!network->Je) { 995883e35e8SHong Zhang ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 99683b2e829SHong Zhang } 997883e35e8SHong Zhang network->Je[3*p] = J[0]; 998883e35e8SHong Zhang network->Je[3*p+1] = J[1]; 999883e35e8SHong Zhang network->Je[3*p+2] = J[2]; 100083b2e829SHong Zhang PetscFunctionReturn(0); 100183b2e829SHong Zhang } 100283b2e829SHong Zhang 100383b2e829SHong Zhang /*@ 100476ddfea5SHong Zhang DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 10051ad426b7SHong Zhang 10061ad426b7SHong Zhang Not Collective 10071ad426b7SHong Zhang 10081ad426b7SHong Zhang Input Parameters: 10091ad426b7SHong Zhang + dm - The DMNetwork object 10101ad426b7SHong Zhang . p - the vertex point 10113e97b6e8SHong Zhang - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 10123e97b6e8SHong Zhang J[0]: this vertex 10133e97b6e8SHong Zhang J[1+2*i]: i-th supporting edge 10143e97b6e8SHong Zhang J[1+2*i+1]: i-th connected vertex 10151ad426b7SHong Zhang 10161ad426b7SHong Zhang Level: intermediate 10171ad426b7SHong Zhang 101883b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix 10191ad426b7SHong Zhang @*/ 1020883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 10215f2c45f1SShri Abhyankar { 10225f2c45f1SShri Abhyankar PetscErrorCode ierr; 10235f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 1024f4431b8cSHong Zhang PetscInt i,*vptr,nedges,vStart,vEnd; 1025883e35e8SHong Zhang const PetscInt *edges; 10265f2c45f1SShri Abhyankar 10275f2c45f1SShri Abhyankar PetscFunctionBegin; 1028883e35e8SHong Zhang if (!network->userVertexJacobian) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 1029883e35e8SHong Zhang 1030883e35e8SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1031f4431b8cSHong Zhang 103283b2e829SHong Zhang if (!network->Jv) { 1033f4431b8cSHong Zhang PetscInt nNodes = network->nNodes,nedges_total; 1034883e35e8SHong Zhang /* count nvertex_total */ 10353e97b6e8SHong Zhang nedges_total = 0; 1036883e35e8SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1037f4431b8cSHong Zhang ierr = PetscMalloc1(nNodes+1,&vptr);CHKERRQ(ierr); 1038f4431b8cSHong Zhang 1039883e35e8SHong Zhang vptr[0] = 0; 1040f4431b8cSHong Zhang for (i=0; i<nNodes; i++) { 1041f4431b8cSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 1042883e35e8SHong Zhang nedges_total += nedges; 1043f4431b8cSHong Zhang vptr[i+1] = vptr[i] + 2*nedges + 1; 10441ad426b7SHong Zhang } 10453e97b6e8SHong Zhang 1046f4431b8cSHong Zhang ierr = PetscCalloc1(2*nedges_total+nNodes,&network->Jv);CHKERRQ(ierr); 1047883e35e8SHong Zhang network->Jvptr = vptr; 1048883e35e8SHong Zhang } 1049883e35e8SHong Zhang 1050883e35e8SHong Zhang vptr = network->Jvptr; 10513e97b6e8SHong Zhang network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 10523e97b6e8SHong Zhang 10533e97b6e8SHong Zhang /* Set Jacobian for each supporting edge and connected vertex */ 1054883e35e8SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 1055883e35e8SHong Zhang for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 1056883e35e8SHong Zhang PetscFunctionReturn(0); 1057883e35e8SHong Zhang } 1058883e35e8SHong Zhang 1059e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 10605cf7da58SHong Zhang { 10615cf7da58SHong Zhang PetscErrorCode ierr; 10625cf7da58SHong Zhang PetscInt j; 10635cf7da58SHong Zhang PetscScalar val=(PetscScalar)ncols; 10645cf7da58SHong Zhang 10655cf7da58SHong Zhang PetscFunctionBegin; 10665cf7da58SHong Zhang if (!ghost) { 10675cf7da58SHong Zhang for (j=0; j<nrows; j++) { 10685cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 10695cf7da58SHong Zhang } 10705cf7da58SHong Zhang } else { 10715cf7da58SHong Zhang for (j=0; j<nrows; j++) { 10725cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 10735cf7da58SHong Zhang } 10745cf7da58SHong Zhang } 10755cf7da58SHong Zhang PetscFunctionReturn(0); 10765cf7da58SHong Zhang } 10775cf7da58SHong Zhang 1078e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 10795cf7da58SHong Zhang { 10805cf7da58SHong Zhang PetscErrorCode ierr; 10815cf7da58SHong Zhang PetscInt j,ncols_u; 10825cf7da58SHong Zhang PetscScalar val; 10835cf7da58SHong Zhang 10845cf7da58SHong Zhang PetscFunctionBegin; 10855cf7da58SHong Zhang if (!ghost) { 10865cf7da58SHong Zhang for (j=0; j<nrows; j++) { 10875cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 10885cf7da58SHong Zhang val = (PetscScalar)ncols_u; 10895cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 10905cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 10915cf7da58SHong Zhang } 10925cf7da58SHong Zhang } else { 10935cf7da58SHong Zhang for (j=0; j<nrows; j++) { 10945cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 10955cf7da58SHong Zhang val = (PetscScalar)ncols_u; 10965cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 10975cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 10985cf7da58SHong Zhang } 10995cf7da58SHong Zhang } 11005cf7da58SHong Zhang PetscFunctionReturn(0); 11015cf7da58SHong Zhang } 11025cf7da58SHong Zhang 1103e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 11045cf7da58SHong Zhang { 11055cf7da58SHong Zhang PetscErrorCode ierr; 11065cf7da58SHong Zhang 11075cf7da58SHong Zhang PetscFunctionBegin; 11085cf7da58SHong Zhang if (Ju) { 11095cf7da58SHong Zhang ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 11105cf7da58SHong Zhang } else { 11115cf7da58SHong Zhang ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 11125cf7da58SHong Zhang } 11135cf7da58SHong Zhang PetscFunctionReturn(0); 11145cf7da58SHong Zhang } 11155cf7da58SHong Zhang 1116e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1117883e35e8SHong Zhang { 1118883e35e8SHong Zhang PetscErrorCode ierr; 1119883e35e8SHong Zhang PetscInt j,*cols; 1120883e35e8SHong Zhang PetscScalar *zeros; 1121883e35e8SHong Zhang 1122883e35e8SHong Zhang PetscFunctionBegin; 1123883e35e8SHong Zhang ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 1124883e35e8SHong Zhang for (j=0; j<ncols; j++) cols[j] = j+ cstart; 1125883e35e8SHong Zhang ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 1126883e35e8SHong Zhang ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 11271ad426b7SHong Zhang PetscFunctionReturn(0); 11281ad426b7SHong Zhang } 1129a4e85ca8SHong Zhang 1130e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 11313e97b6e8SHong Zhang { 11323e97b6e8SHong Zhang PetscErrorCode ierr; 11333e97b6e8SHong Zhang PetscInt j,M,N,row,col,ncols_u; 11343e97b6e8SHong Zhang const PetscInt *cols; 11353e97b6e8SHong Zhang PetscScalar zero=0.0; 11363e97b6e8SHong Zhang 11373e97b6e8SHong Zhang PetscFunctionBegin; 11383e97b6e8SHong Zhang ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 11393e97b6e8SHong 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); 11403e97b6e8SHong Zhang 11413e97b6e8SHong Zhang for (row=0; row<nrows; row++) { 11423e97b6e8SHong Zhang ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 11433e97b6e8SHong Zhang for (j=0; j<ncols_u; j++) { 11443e97b6e8SHong Zhang col = cols[j] + cstart; 11453e97b6e8SHong Zhang ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 11463e97b6e8SHong Zhang } 11473e97b6e8SHong Zhang ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 11483e97b6e8SHong Zhang } 11493e97b6e8SHong Zhang PetscFunctionReturn(0); 11503e97b6e8SHong Zhang } 11511ad426b7SHong Zhang 1152e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1153a4e85ca8SHong Zhang { 1154a4e85ca8SHong Zhang PetscErrorCode ierr; 1155f4431b8cSHong Zhang 1156a4e85ca8SHong Zhang PetscFunctionBegin; 1157a4e85ca8SHong Zhang if (Ju) { 1158a4e85ca8SHong Zhang ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1159a4e85ca8SHong Zhang } else { 1160a4e85ca8SHong Zhang ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1161a4e85ca8SHong Zhang } 1162a4e85ca8SHong Zhang PetscFunctionReturn(0); 1163a4e85ca8SHong Zhang } 1164a4e85ca8SHong Zhang 116524121865SAdrian Maldonado 116624121865SAdrian 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. 116724121865SAdrian Maldonado */ 116824121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 116924121865SAdrian Maldonado { 117024121865SAdrian Maldonado PetscErrorCode ierr; 117124121865SAdrian Maldonado PetscInt i, size, dof; 117224121865SAdrian Maldonado PetscInt *glob2loc; 117324121865SAdrian Maldonado 117424121865SAdrian Maldonado PetscFunctionBegin; 117524121865SAdrian Maldonado ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr); 117624121865SAdrian Maldonado ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr); 117724121865SAdrian Maldonado 117824121865SAdrian Maldonado for (i = 0; i < size; i++) { 117924121865SAdrian Maldonado ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr); 118024121865SAdrian Maldonado dof = (dof >= 0) ? dof : -(dof + 1); 118124121865SAdrian Maldonado glob2loc[i] = dof; 118224121865SAdrian Maldonado } 118324121865SAdrian Maldonado 118424121865SAdrian Maldonado ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 118524121865SAdrian Maldonado #if 0 118624121865SAdrian Maldonado ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 118724121865SAdrian Maldonado #endif 118824121865SAdrian Maldonado PetscFunctionReturn(0); 118924121865SAdrian Maldonado } 119024121865SAdrian Maldonado 11911ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 11921ad426b7SHong Zhang { 11931ad426b7SHong Zhang PetscErrorCode ierr; 119424121865SAdrian Maldonado PetscMPIInt rank, size; 11951ad426b7SHong Zhang DM_Network *network = (DM_Network*) dm->data; 1196a4e85ca8SHong Zhang PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 1197840c2264SHong Zhang PetscInt cstart,ncols,j,e,v; 119824121865SAdrian Maldonado PetscBool ghost,ghost_vc,ghost2,isNest; 1199a4e85ca8SHong Zhang Mat Juser; 1200bfbc38dcSHong Zhang PetscSection sectionGlobal; 1201447d78afSSatish Balay PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 1202a4e85ca8SHong Zhang const PetscInt *edges,*cone; 12035cf7da58SHong Zhang MPI_Comm comm; 120424121865SAdrian Maldonado MatType mtype; 12055cf7da58SHong Zhang Vec vd_nz,vo_nz; 12065cf7da58SHong Zhang PetscInt *dnnz,*onnz; 12075cf7da58SHong Zhang PetscScalar *vdnz,*vonz; 12081ad426b7SHong Zhang 12091ad426b7SHong Zhang PetscFunctionBegin; 121024121865SAdrian Maldonado mtype = dm->mattype; 121124121865SAdrian Maldonado ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr); 121224121865SAdrian Maldonado 121324121865SAdrian Maldonado if (isNest) { 12140731d606SHong Zhang /* ierr = DMCreateMatrix_Network_Nest(); */ 121524121865SAdrian Maldonado PetscInt eDof, vDof; 121624121865SAdrian Maldonado Mat j11, j12, j21, j22, bA[2][2]; 121724121865SAdrian Maldonado ISLocalToGlobalMapping eISMap, vISMap; 121824121865SAdrian Maldonado 121924121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 122024121865SAdrian Maldonado ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 122124121865SAdrian Maldonado ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 122224121865SAdrian Maldonado 122324121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr); 122424121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr); 122524121865SAdrian Maldonado 122624121865SAdrian Maldonado ierr = MatCreate(PETSC_COMM_WORLD, &j11);CHKERRQ(ierr); 122724121865SAdrian Maldonado ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 122824121865SAdrian Maldonado ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr); 122924121865SAdrian Maldonado 123024121865SAdrian Maldonado ierr = MatCreate(PETSC_COMM_WORLD, &j12);CHKERRQ(ierr); 123124121865SAdrian Maldonado ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr); 123224121865SAdrian Maldonado ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr); 123324121865SAdrian Maldonado 123424121865SAdrian Maldonado ierr = MatCreate(PETSC_COMM_WORLD, &j21);CHKERRQ(ierr); 123524121865SAdrian Maldonado ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 123624121865SAdrian Maldonado ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr); 123724121865SAdrian Maldonado 123824121865SAdrian Maldonado ierr = MatCreate(PETSC_COMM_WORLD, &j22);CHKERRQ(ierr); 123924121865SAdrian Maldonado ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 124024121865SAdrian Maldonado ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr); 124124121865SAdrian Maldonado 124224121865SAdrian Maldonado bA[0][0] = j11; 124324121865SAdrian Maldonado bA[0][1] = j12; 124424121865SAdrian Maldonado bA[1][0] = j21; 124524121865SAdrian Maldonado bA[1][1] = j22; 124624121865SAdrian Maldonado 124724121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr); 124824121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr); 124924121865SAdrian Maldonado 125024121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr); 125124121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr); 125224121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr); 125324121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr); 125424121865SAdrian Maldonado 125524121865SAdrian Maldonado ierr = MatSetUp(j11);CHKERRQ(ierr); 125624121865SAdrian Maldonado ierr = MatSetUp(j12);CHKERRQ(ierr); 125724121865SAdrian Maldonado ierr = MatSetUp(j21);CHKERRQ(ierr); 125824121865SAdrian Maldonado ierr = MatSetUp(j22);CHKERRQ(ierr); 125924121865SAdrian Maldonado 126024121865SAdrian Maldonado ierr = MatCreateNest(PETSC_COMM_WORLD,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr); 126124121865SAdrian Maldonado ierr = MatSetUp(*J);CHKERRQ(ierr); 126224121865SAdrian Maldonado ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr); 126324121865SAdrian Maldonado ierr = MatDestroy(&j11);CHKERRQ(ierr); 126424121865SAdrian Maldonado ierr = MatDestroy(&j12);CHKERRQ(ierr); 126524121865SAdrian Maldonado ierr = MatDestroy(&j21);CHKERRQ(ierr); 126624121865SAdrian Maldonado ierr = MatDestroy(&j22);CHKERRQ(ierr); 126724121865SAdrian Maldonado 126824121865SAdrian Maldonado ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 126924121865SAdrian Maldonado ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 127024121865SAdrian Maldonado 127124121865SAdrian Maldonado /* Free structures */ 127224121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr); 127324121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr); 127424121865SAdrian Maldonado 127524121865SAdrian Maldonado PetscFunctionReturn(0); 127624121865SAdrian Maldonado } else if (!network->userEdgeJacobian && !network->userVertexJacobian) { 1277a4e85ca8SHong Zhang /* user does not provide Jacobian blocks */ 1278bfbc38dcSHong Zhang ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 1279bfbc38dcSHong Zhang ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 12801ad426b7SHong Zhang PetscFunctionReturn(0); 12811ad426b7SHong Zhang } 12821ad426b7SHong Zhang 1283bfbc38dcSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 12842a945128SHong Zhang ierr = DMGetDefaultGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 1285bfbc38dcSHong Zhang ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 1286bfbc38dcSHong Zhang ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 12872a945128SHong Zhang 12882a945128SHong Zhang ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 12892a945128SHong Zhang ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 129089898e50SHong Zhang 129189898e50SHong Zhang /* (1) Set matrix preallocation */ 129289898e50SHong Zhang /*------------------------------*/ 1293840c2264SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1294840c2264SHong Zhang ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr); 1295840c2264SHong Zhang ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr); 1296840c2264SHong Zhang ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr); 1297840c2264SHong Zhang ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr); 1298840c2264SHong Zhang ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr); 1299840c2264SHong Zhang 130089898e50SHong Zhang /* Set preallocation for edges */ 130189898e50SHong Zhang /*-----------------------------*/ 1302840c2264SHong Zhang ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 1303840c2264SHong Zhang 1304bdcb62a2SHong Zhang ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr); 1305840c2264SHong Zhang for (e=eStart; e<eEnd; e++) { 1306840c2264SHong Zhang /* Get row indices */ 1307840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 1308840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 1309840c2264SHong Zhang if (nrows) { 1310840c2264SHong Zhang if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je"); 1311840c2264SHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 1312840c2264SHong Zhang 13135cf7da58SHong Zhang /* Set preallocation for conntected vertices */ 1314840c2264SHong Zhang ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr); 1315840c2264SHong Zhang for (v=0; v<2; v++) { 1316840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 1317840c2264SHong Zhang 1318840c2264SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 1319840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 13205cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1321840c2264SHong Zhang } 1322840c2264SHong Zhang 132389898e50SHong Zhang /* Set preallocation for edge self */ 1324840c2264SHong Zhang cstart = rstart; 1325840c2264SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 13265cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1327840c2264SHong Zhang } 1328840c2264SHong Zhang } 1329840c2264SHong Zhang 133089898e50SHong Zhang /* Set preallocation for vertices */ 133189898e50SHong Zhang /*--------------------------------*/ 1332840c2264SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 1333840c2264SHong Zhang if (vEnd - vStart) { 1334840c2264SHong Zhang if (!network->Jv) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Jv"); 1335840c2264SHong Zhang vptr = network->Jvptr; 1336840c2264SHong Zhang if (!vptr) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide vptr"); 1337840c2264SHong Zhang } 1338840c2264SHong Zhang 1339840c2264SHong Zhang for (v=vStart; v<vEnd; v++) { 1340840c2264SHong Zhang /* Get row indices */ 1341840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1342840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 1343840c2264SHong Zhang if (!nrows) continue; 1344840c2264SHong Zhang 1345bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1346bdcb62a2SHong Zhang if (ghost) { 1347bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1348bdcb62a2SHong Zhang } else { 1349bdcb62a2SHong Zhang rows_v = rows; 1350bdcb62a2SHong Zhang } 1351bdcb62a2SHong Zhang 1352bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1353840c2264SHong Zhang 1354840c2264SHong Zhang /* Get supporting edges and connected vertices */ 1355840c2264SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1356840c2264SHong Zhang 1357840c2264SHong Zhang for (e=0; e<nedges; e++) { 1358840c2264SHong Zhang /* Supporting edges */ 1359840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1360840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1361840c2264SHong Zhang 1362840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1363bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1364840c2264SHong Zhang 1365840c2264SHong Zhang /* Connected vertices */ 1366840c2264SHong Zhang ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr); 1367840c2264SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 1368840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr); 1369840c2264SHong Zhang 1370840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1371840c2264SHong Zhang 1372840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1373e102a522SHong Zhang if (ghost_vc||ghost) { 1374e102a522SHong Zhang ghost2 = PETSC_TRUE; 1375e102a522SHong Zhang } else { 1376e102a522SHong Zhang ghost2 = PETSC_FALSE; 1377e102a522SHong Zhang } 1378e102a522SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr); 1379840c2264SHong Zhang } 1380840c2264SHong Zhang 138189898e50SHong Zhang /* Set preallocation for vertex self */ 1382840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1383840c2264SHong Zhang if (!ghost) { 1384840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1385840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1386bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1387840c2264SHong Zhang } 1388bdcb62a2SHong Zhang if (ghost) { 1389bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1390bdcb62a2SHong Zhang } 1391840c2264SHong Zhang } 1392840c2264SHong Zhang 1393840c2264SHong Zhang ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr); 1394840c2264SHong Zhang ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr); 13955cf7da58SHong Zhang 13965cf7da58SHong Zhang ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr); 13975cf7da58SHong Zhang 13985cf7da58SHong Zhang ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr); 1399840c2264SHong Zhang ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr); 1400840c2264SHong Zhang 1401840c2264SHong Zhang ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr); 1402840c2264SHong Zhang ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr); 1403840c2264SHong Zhang for (j=0; j<localSize; j++) { 1404e102a522SHong Zhang dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 1405e102a522SHong Zhang onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 1406840c2264SHong Zhang } 1407840c2264SHong Zhang ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr); 1408840c2264SHong Zhang ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr); 1409840c2264SHong Zhang ierr = VecDestroy(&vd_nz);CHKERRQ(ierr); 1410840c2264SHong Zhang ierr = VecDestroy(&vo_nz);CHKERRQ(ierr); 1411840c2264SHong Zhang 14125cf7da58SHong Zhang ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr); 14135cf7da58SHong Zhang ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr); 14145cf7da58SHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 14155cf7da58SHong Zhang 14165cf7da58SHong Zhang ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr); 14175cf7da58SHong Zhang 141889898e50SHong Zhang /* (2) Set matrix entries for edges */ 141989898e50SHong Zhang /*----------------------------------*/ 14201ad426b7SHong Zhang for (e=eStart; e<eEnd; e++) { 1421bfbc38dcSHong Zhang /* Get row indices */ 14221ad426b7SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 142317df6e9eSHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 14244b976069SHong Zhang if (nrows) { 14254b976069SHong Zhang if (!network->Je) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Uer must provide Je"); 1426bdcb62a2SHong Zhang 142717df6e9eSHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 14281ad426b7SHong Zhang 1429bfbc38dcSHong Zhang /* Set matrix entries for conntected vertices */ 143017df6e9eSHong Zhang ierr = DMNetworkGetConnectedNodes(dm,e,&cone);CHKERRQ(ierr); 1431bfbc38dcSHong Zhang for (v=0; v<2; v++) { 1432bfbc38dcSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 1433883e35e8SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 14343e97b6e8SHong Zhang 1435a4e85ca8SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 1436a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1437bfbc38dcSHong Zhang } 143817df6e9eSHong Zhang 1439bfbc38dcSHong Zhang /* Set matrix entries for edge self */ 14403e97b6e8SHong Zhang cstart = rstart; 1441a4e85ca8SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 1442a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 14431ad426b7SHong Zhang } 14444b976069SHong Zhang } 14451ad426b7SHong Zhang 1446bfbc38dcSHong Zhang /* Set matrix entries for vertices */ 144783b2e829SHong Zhang /*---------------------------------*/ 14481ad426b7SHong Zhang for (v=vStart; v<vEnd; v++) { 1449bfbc38dcSHong Zhang /* Get row indices */ 1450596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1451596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 14524b976069SHong Zhang if (!nrows) continue; 1453596e729fSHong Zhang 1454bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1455bdcb62a2SHong Zhang if (ghost) { 1456bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1457bdcb62a2SHong Zhang } else { 1458bdcb62a2SHong Zhang rows_v = rows; 1459bdcb62a2SHong Zhang } 1460bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1461596e729fSHong Zhang 1462bfbc38dcSHong Zhang /* Get supporting edges and connected vertices */ 1463596e729fSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1464596e729fSHong Zhang 1465596e729fSHong Zhang for (e=0; e<nedges; e++) { 1466bfbc38dcSHong Zhang /* Supporting edges */ 1467596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1468596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1469596e729fSHong Zhang 1470a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 1471bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1472596e729fSHong Zhang 1473bfbc38dcSHong Zhang /* Connected vertices */ 147444aca652SHong Zhang ierr = DMNetworkGetConnectedNodes(dm,edges[e],&cone);CHKERRQ(ierr); 14752a945128SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 14762a945128SHong Zhang 147744aca652SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 147844aca652SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1479a4e85ca8SHong Zhang 1480a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 1481bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1482596e729fSHong Zhang } 1483596e729fSHong Zhang 1484bfbc38dcSHong Zhang /* Set matrix entries for vertex self */ 14851ad426b7SHong Zhang if (!ghost) { 1486596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 1487a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 1488bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr); 1489bdcb62a2SHong Zhang } 1490bdcb62a2SHong Zhang if (ghost) { 1491bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1492bdcb62a2SHong Zhang } 14931ad426b7SHong Zhang } 1494a4e85ca8SHong Zhang ierr = PetscFree(rows);CHKERRQ(ierr); 1495bdcb62a2SHong Zhang 14961ad426b7SHong Zhang ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 14971ad426b7SHong Zhang ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1498dd6f46cdSHong Zhang 14995f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 15005f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15015f2c45f1SShri Abhyankar } 15025f2c45f1SShri Abhyankar 15035f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 15045f2c45f1SShri Abhyankar { 15055f2c45f1SShri Abhyankar PetscErrorCode ierr; 15065f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15075f2c45f1SShri Abhyankar 15085f2c45f1SShri Abhyankar PetscFunctionBegin; 15098415c774SShri Abhyankar if (--network->refct > 0) PetscFunctionReturn(0); 151083b2e829SHong Zhang if (network->Je) { 151183b2e829SHong Zhang ierr = PetscFree(network->Je);CHKERRQ(ierr); 151283b2e829SHong Zhang } 151383b2e829SHong Zhang if (network->Jv) { 1514883e35e8SHong Zhang ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 151583b2e829SHong Zhang ierr = PetscFree(network->Jv);CHKERRQ(ierr); 15161ad426b7SHong Zhang } 151713c2a604SAdrian Maldonado 151813c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr); 151913c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr); 152013c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr); 152113c2a604SAdrian Maldonado if (network->vertex.sf) { 152213c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr); 152313c2a604SAdrian Maldonado } 152413c2a604SAdrian Maldonado /* edge */ 152513c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr); 152613c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr); 152713c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr); 152813c2a604SAdrian Maldonado if (network->edge.sf) { 152913c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr); 153013c2a604SAdrian Maldonado } 15315f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 15325f2c45f1SShri Abhyankar network->edges = NULL; 15335f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 15345f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 153583b2e829SHong Zhang 15365f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 15375f2c45f1SShri Abhyankar ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 15385f2c45f1SShri Abhyankar ierr = PetscFree(network->header);CHKERRQ(ierr); 15395f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 15405f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15415f2c45f1SShri Abhyankar } 15425f2c45f1SShri Abhyankar 15435f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 15445f2c45f1SShri Abhyankar { 15455f2c45f1SShri Abhyankar PetscErrorCode ierr; 15465f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15475f2c45f1SShri Abhyankar 15485f2c45f1SShri Abhyankar PetscFunctionBegin; 15495f2c45f1SShri Abhyankar ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 15505f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15515f2c45f1SShri Abhyankar } 15525f2c45f1SShri Abhyankar 15535f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 15545f2c45f1SShri Abhyankar { 15555f2c45f1SShri Abhyankar PetscErrorCode ierr; 15565f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15575f2c45f1SShri Abhyankar 15585f2c45f1SShri Abhyankar PetscFunctionBegin; 15595f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 15605f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15615f2c45f1SShri Abhyankar } 15625f2c45f1SShri Abhyankar 15635f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 15645f2c45f1SShri Abhyankar { 15655f2c45f1SShri Abhyankar PetscErrorCode ierr; 15665f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15675f2c45f1SShri Abhyankar 15685f2c45f1SShri Abhyankar PetscFunctionBegin; 15695f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 15705f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15715f2c45f1SShri Abhyankar } 15725f2c45f1SShri Abhyankar 15735f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 15745f2c45f1SShri Abhyankar { 15755f2c45f1SShri Abhyankar PetscErrorCode ierr; 15765f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15775f2c45f1SShri Abhyankar 15785f2c45f1SShri Abhyankar PetscFunctionBegin; 15795f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 15805f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15815f2c45f1SShri Abhyankar } 15825f2c45f1SShri Abhyankar 15835f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 15845f2c45f1SShri Abhyankar { 15855f2c45f1SShri Abhyankar PetscErrorCode ierr; 15865f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 15875f2c45f1SShri Abhyankar 15885f2c45f1SShri Abhyankar PetscFunctionBegin; 15895f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 15905f2c45f1SShri Abhyankar PetscFunctionReturn(0); 15915f2c45f1SShri Abhyankar } 1592