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 /*@ 28e2aaf10cSShri Abhyankar DMNetworkSetSizes - Sets the number of subnetworks,local and global vertices and edges for each subnetwork. 295f2c45f1SShri Abhyankar 305f2c45f1SShri Abhyankar Collective on DM 315f2c45f1SShri Abhyankar 325f2c45f1SShri Abhyankar Input Parameters: 335f2c45f1SShri Abhyankar + dm - the dm object 34e2aaf10cSShri Abhyankar . Nsubnet - number of subnetworks 35e2aaf10cSShri Abhyankar . nV - number of local vertices for each subnetwork 36e2aaf10cSShri Abhyankar . nE - number of local edges for each subnetwork 37e2aaf10cSShri Abhyankar . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork 38e2aaf10cSShri Abhyankar - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork 395f2c45f1SShri Abhyankar 405f2c45f1SShri Abhyankar Notes 415f2c45f1SShri Abhyankar If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 425f2c45f1SShri Abhyankar 435f2c45f1SShri Abhyankar You cannot change the sizes once they have been set 445f2c45f1SShri Abhyankar 451b266c99SBarry Smith Level: intermediate 461b266c99SBarry Smith 471b266c99SBarry Smith .seealso: DMNetworkCreate() 485f2c45f1SShri Abhyankar @*/ 49e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetSizes(DM dm, PetscInt Nsubnet, PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[]) 505f2c45f1SShri Abhyankar { 515f2c45f1SShri Abhyankar PetscErrorCode ierr; 525f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 53e2aaf10cSShri Abhyankar PetscInt a[2],b[2],i; 545f2c45f1SShri Abhyankar 555f2c45f1SShri Abhyankar PetscFunctionBegin; 565f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 57e2aaf10cSShri Abhyankar if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet); 58e2aaf10cSShri Abhyankar 59e2aaf10cSShri Abhyankar if(Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2); 60e2aaf10cSShri Abhyankar if(network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network"); 61e2aaf10cSShri Abhyankar 62e2aaf10cSShri Abhyankar network->nsubnet = Nsubnet; 63e2aaf10cSShri Abhyankar ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr); 64e2aaf10cSShri Abhyankar for(i=0; i < network->nsubnet; i++) { 65e2aaf10cSShri Abhyankar if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],5); 66e2aaf10cSShri Abhyankar if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],6); 67e2aaf10cSShri Abhyankar if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]); 68e2aaf10cSShri Abhyankar if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]); 69e2aaf10cSShri Abhyankar a[0] = nV[i]; a[1] = nE[i]; 70b2566f29SBarry Smith ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 71e2aaf10cSShri Abhyankar network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1]; 72e2aaf10cSShri Abhyankar 73e2aaf10cSShri Abhyankar network->subnet[i].id = i; 74e2aaf10cSShri Abhyankar 75e2aaf10cSShri Abhyankar network->subnet[i].nvtx = nV[i]; 76e2aaf10cSShri Abhyankar network->subnet[i].vStart = network->nVertices; 77e2aaf10cSShri Abhyankar network->subnet[i].vEnd = network->subnet[i].vStart + network->subnet[i].Nvtx; 78e2aaf10cSShri Abhyankar network->nVertices += network->subnet[i].nvtx; 79e2aaf10cSShri Abhyankar network->NVertices += network->subnet[i].Nvtx; 80e2aaf10cSShri Abhyankar 81e2aaf10cSShri Abhyankar network->subnet[i].nedge = nE[i]; 82e2aaf10cSShri Abhyankar network->subnet[i].eStart = network->nEdges; 83e2aaf10cSShri Abhyankar network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge; 84e2aaf10cSShri Abhyankar network->nEdges += network->subnet[i].nedge; 85e2aaf10cSShri Abhyankar network->NEdges += network->subnet[i].Nedge; 865f2c45f1SShri Abhyankar } 875f2c45f1SShri Abhyankar PetscFunctionReturn(0); 885f2c45f1SShri Abhyankar } 895f2c45f1SShri Abhyankar 907765340cSHong Zhang PetscErrorCode DMNetworkSetSizesCoupled(DM dm,PetscInt Nsubnet,PetscInt NsubnetCouple,PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[]) 917765340cSHong Zhang { 927765340cSHong Zhang PetscErrorCode ierr; 937765340cSHong Zhang DM_Network *network = (DM_Network*) dm->data; 947765340cSHong Zhang PetscInt a[2],b[2],i; 957765340cSHong Zhang 967765340cSHong Zhang PetscFunctionBegin; 97*34bd1bf5SHong Zhang //printf("DMNetworkSetSizesCoupled...Nsubnet %d, NsubnetCouple %d\n",Nsubnet,NsubnetCouple); 987765340cSHong Zhang PetscValidHeaderSpecific(dm,DM_CLASSID,1); 997765340cSHong Zhang if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet); 1007765340cSHong Zhang if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple); 1017765340cSHong Zhang 1027765340cSHong Zhang if (Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2); 1037765340cSHong Zhang if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network"); 1047765340cSHong Zhang 1057765340cSHong Zhang network->nsubnet = Nsubnet + NsubnetCouple; 1067765340cSHong Zhang Nsubnet += NsubnetCouple; 1077765340cSHong Zhang network->ncsubnet = NsubnetCouple; 1087765340cSHong Zhang ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr); 1097765340cSHong Zhang for(i=0; i < Nsubnet; i++) { 1107765340cSHong Zhang if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],6); 1117765340cSHong Zhang if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],7); 1127765340cSHong Zhang if (NV[i] > 0 && nV[i] > NV[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local vertex size %D cannot be larger than global vertex size %D",i,nV[i],NV[i]); 1137765340cSHong Zhang if (NE[i] > 0 && nE[i] > NE[i]) SETERRQ3(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Subnetwork %D: Local edge size %D cannot be larger than global edge size %D",i,nE[i],NE[i]); 1147765340cSHong Zhang a[0] = nV[i]; a[1] = nE[i]; 1157765340cSHong Zhang ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 1167765340cSHong Zhang network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1]; 1177765340cSHong Zhang 1187765340cSHong Zhang network->subnet[i].id = i; 1197765340cSHong Zhang 1207765340cSHong Zhang network->subnet[i].nvtx = nV[i]; 1217765340cSHong Zhang network->subnet[i].vStart = network->nVertices; 1227765340cSHong Zhang network->subnet[i].vEnd = network->subnet[i].vStart + network->subnet[i].Nvtx; 1237765340cSHong Zhang network->nVertices += network->subnet[i].nvtx; 1247765340cSHong Zhang network->NVertices += network->subnet[i].Nvtx; 1257765340cSHong Zhang 1267765340cSHong Zhang network->subnet[i].nedge = nE[i]; 1277765340cSHong Zhang network->subnet[i].eStart = network->nEdges; 1287765340cSHong Zhang network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge; 1297765340cSHong Zhang network->nEdges += network->subnet[i].nedge; 1307765340cSHong Zhang network->NEdges += network->subnet[i].Nedge; 1317765340cSHong Zhang } 1327765340cSHong Zhang PetscFunctionReturn(0); 1337765340cSHong Zhang } 1347765340cSHong Zhang 1355f2c45f1SShri Abhyankar /*@ 1365f2c45f1SShri Abhyankar DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 1375f2c45f1SShri Abhyankar 1385f2c45f1SShri Abhyankar Logically collective on DM 1395f2c45f1SShri Abhyankar 1405f2c45f1SShri Abhyankar Input Parameters: 141e2aaf10cSShri Abhyankar . edges - list of edges for each subnetwork 1425f2c45f1SShri Abhyankar 1435f2c45f1SShri Abhyankar Notes: 1445f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 1455f2c45f1SShri Abhyankar not be destroyed before the call to DMNetworkLayoutSetUp 1465f2c45f1SShri Abhyankar 1475f2c45f1SShri Abhyankar Level: intermediate 1485f2c45f1SShri Abhyankar 1495f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes 1505f2c45f1SShri Abhyankar @*/ 151e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int *edgelist[]) 1525f2c45f1SShri Abhyankar { 1535f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 154e2aaf10cSShri Abhyankar PetscInt i; 1555f2c45f1SShri Abhyankar 1565f2c45f1SShri Abhyankar PetscFunctionBegin; 157e2aaf10cSShri Abhyankar for(i=0; i < network->nsubnet; i++) { 158e2aaf10cSShri Abhyankar network->subnet[i].edgelist = edgelist[i]; 159e2aaf10cSShri Abhyankar } 1605f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1615f2c45f1SShri Abhyankar } 1625f2c45f1SShri Abhyankar 1635f2c45f1SShri Abhyankar /*@ 1645f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 1655f2c45f1SShri Abhyankar 1665f2c45f1SShri Abhyankar Collective on DM 1675f2c45f1SShri Abhyankar 1685f2c45f1SShri Abhyankar Input Parameters 1695f2c45f1SShri Abhyankar . DM - the dmnetwork object 1705f2c45f1SShri Abhyankar 1715f2c45f1SShri Abhyankar Notes: 1725f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 1735f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 1745f2c45f1SShri Abhyankar 1755f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 1765f2c45f1SShri Abhyankar 1775f2c45f1SShri Abhyankar Level: intermediate 1785f2c45f1SShri Abhyankar 1795f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 1805f2c45f1SShri Abhyankar @*/ 1815f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 1825f2c45f1SShri Abhyankar { 1835f2c45f1SShri Abhyankar PetscErrorCode ierr; 1845f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1855f2c45f1SShri Abhyankar PetscInt dim = 1; /* One dimensional network */ 1865f2c45f1SShri Abhyankar PetscInt numCorners=2; 1875f2c45f1SShri Abhyankar PetscInt spacedim=2; 1885f2c45f1SShri Abhyankar double *vertexcoords=NULL; 189e2aaf10cSShri Abhyankar PetscInt i,j; 1905f2c45f1SShri Abhyankar PetscInt ndata; 191e2aaf10cSShri Abhyankar PetscInt ctr=0; 1925f2c45f1SShri Abhyankar 1935f2c45f1SShri Abhyankar PetscFunctionBegin; 1946fefedf4SHong Zhang if (network->nVertices) { 1956fefedf4SHong Zhang ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr); 1965f2c45f1SShri Abhyankar } 197e2aaf10cSShri Abhyankar 198e2aaf10cSShri Abhyankar /* Create the edgelist for the network by concatenating edgelists of the subnetworks */ 199e2aaf10cSShri Abhyankar ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr); 200e2aaf10cSShri Abhyankar for(i=0; i < network->nsubnet; i++) { 201e2aaf10cSShri Abhyankar for(j = 0; j < network->subnet[i].nedge; j++) { 202e2aaf10cSShri Abhyankar network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j]; 203e2aaf10cSShri Abhyankar network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1]; 204e2aaf10cSShri Abhyankar ctr++; 205e2aaf10cSShri Abhyankar } 206e2aaf10cSShri Abhyankar } 207e2aaf10cSShri Abhyankar 20861de3474SHong Zhang #if 0 209e2aaf10cSShri Abhyankar for(i=0; i < network->nEdges; i++) { 210e2aaf10cSShri Abhyankar ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr); 211e2aaf10cSShri Abhyankar } 212e2aaf10cSShri Abhyankar #endif 213e2aaf10cSShri Abhyankar 2146fefedf4SHong Zhang ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 2156fefedf4SHong Zhang if (network->nVertices) { 2165f2c45f1SShri Abhyankar ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 2175f2c45f1SShri Abhyankar } 218e2aaf10cSShri Abhyankar ierr = PetscFree(network->edges);CHKERRQ(ierr); 219e2aaf10cSShri Abhyankar 2205f2c45f1SShri Abhyankar ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 2215f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 2225f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 2235f2c45f1SShri Abhyankar 2245f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 2255f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 2265f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 2275f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 2285f2c45f1SShri Abhyankar 2292727e31bSShri Abhyankar /* Create vertices and edges array for the subnetworks */ 2302727e31bSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 2312727e31bSShri Abhyankar ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr); 2322727e31bSShri Abhyankar ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr); 2332727e31bSShri Abhyankar /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 2342727e31bSShri Abhyankar These get updated when the vertices and edges are added. */ 2352727e31bSShri Abhyankar network->subnet[j].nvtx = network->subnet[j].nedge = 0; 2362727e31bSShri Abhyankar } 2372727e31bSShri Abhyankar 2385f2c45f1SShri Abhyankar network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 2396caa05f4SBarry Smith ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 240e2aaf10cSShri Abhyankar for(i=network->eStart; i < network->eEnd; i++) { 241e2aaf10cSShri Abhyankar network->header[i].index = i; /* Global edge number */ 242e2aaf10cSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 243e2aaf10cSShri Abhyankar if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) { 244e2aaf10cSShri Abhyankar network->header[i].subnetid = j; /* Subnetwork id */ 2452727e31bSShri Abhyankar network->subnet[j].edges[network->subnet[j].nedge++] = i; 246e2aaf10cSShri Abhyankar break; 247e2aaf10cSShri Abhyankar } 2487b6afd5bSHong Zhang } 2495f2c45f1SShri Abhyankar network->header[i].ndata = 0; 2505f2c45f1SShri Abhyankar ndata = network->header[i].ndata; 2515f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 2525f2c45f1SShri Abhyankar network->header[i].offset[ndata] = 0; 2535f2c45f1SShri Abhyankar } 254e2aaf10cSShri Abhyankar 255e2aaf10cSShri Abhyankar for(i=network->vStart; i < network->vEnd; i++) { 256e2aaf10cSShri Abhyankar network->header[i].index = i - network->vStart; 257e2aaf10cSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 258e2aaf10cSShri Abhyankar if((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) { 259e2aaf10cSShri Abhyankar network->header[i].subnetid = j; 2602727e31bSShri Abhyankar network->subnet[j].vertices[network->subnet[j].nvtx++] = i; 261e2aaf10cSShri Abhyankar break; 262e2aaf10cSShri Abhyankar } 263e2aaf10cSShri Abhyankar } 264e2aaf10cSShri Abhyankar network->header[i].ndata = 0; 265e2aaf10cSShri Abhyankar ndata = network->header[i].ndata; 266e2aaf10cSShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 267e2aaf10cSShri Abhyankar network->header[i].offset[ndata] = 0; 268e2aaf10cSShri Abhyankar } 269e2aaf10cSShri Abhyankar 270854ce69bSBarry Smith ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 2716500d4abSHong Zhang PetscFunctionReturn(0); 2726500d4abSHong Zhang } 273e2aaf10cSShri Abhyankar 2746500d4abSHong Zhang PetscErrorCode DMNetworkLayoutSetUpCoupled(DM dm) 2756500d4abSHong Zhang { 2766500d4abSHong Zhang PetscErrorCode ierr; 2776500d4abSHong Zhang DM_Network *network = (DM_Network*)dm->data; 2786500d4abSHong Zhang PetscInt dim = 1; /* One dimensional network */ 279991cf414SHong Zhang PetscInt numCorners=2,spacedim=2; 2806500d4abSHong Zhang double *vertexcoords=NULL; 2817765340cSHong Zhang PetscInt i,j,ndata,ctr=0,nsubnet; 282991cf414SHong Zhang PetscInt *edgelist_couple=NULL,k,netid,vid; 2836500d4abSHong Zhang 2846500d4abSHong Zhang PetscFunctionBegin; 2856500d4abSHong Zhang if (network->nVertices) { 2866500d4abSHong Zhang ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr); 2876500d4abSHong Zhang } 2886500d4abSHong Zhang 2896500d4abSHong Zhang /* Create the edgelist for the network by concatenating edgelists of the subnetworks */ 2907765340cSHong Zhang nsubnet = network->nsubnet - network->ncsubnet; 2916500d4abSHong Zhang ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr); 2927765340cSHong Zhang for (i=0; i < nsubnet; i++) { 2936500d4abSHong Zhang for (j = 0; j < network->subnet[i].nedge; j++) { 2946500d4abSHong Zhang network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j]; 2956500d4abSHong Zhang network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1]; 2966500d4abSHong Zhang ctr++; 2976500d4abSHong Zhang } 2986500d4abSHong Zhang } 2997765340cSHong Zhang 3007765340cSHong Zhang i = nsubnet; /* netid of coupling subnet */ 3017765340cSHong Zhang nsubnet = network->nsubnet; 3027765340cSHong Zhang while (i < nsubnet) { 303991cf414SHong Zhang edgelist_couple = network->subnet[i].edgelist; 3046500d4abSHong Zhang k = 0; 3056500d4abSHong Zhang for (j = 0; j < network->subnet[i].nedge; j++) { 3066500d4abSHong Zhang netid = edgelist_couple[k]; vid = edgelist_couple[k+1]; 3076500d4abSHong Zhang network->edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2; 3086500d4abSHong Zhang 3096500d4abSHong Zhang netid = edgelist_couple[k]; vid = edgelist_couple[k+1]; 310991cf414SHong Zhang network->edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2; 3116500d4abSHong Zhang ctr++; 3126500d4abSHong Zhang } 3137765340cSHong Zhang i++; 3147765340cSHong Zhang } 3156500d4abSHong Zhang 316*34bd1bf5SHong Zhang #if 0 3176500d4abSHong Zhang for(i=0; i < network->nEdges; i++) { 3186500d4abSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr); 3196500d4abSHong Zhang printf("\n"); 3206500d4abSHong Zhang } 3216500d4abSHong Zhang #endif 3226500d4abSHong Zhang 3236500d4abSHong Zhang ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 3246500d4abSHong Zhang if (network->nVertices) { 3256500d4abSHong Zhang ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 3266500d4abSHong Zhang } 3276500d4abSHong Zhang ierr = PetscFree(network->edges);CHKERRQ(ierr); 3286500d4abSHong Zhang 3296500d4abSHong Zhang ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 3306500d4abSHong Zhang ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 3316500d4abSHong Zhang ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 3326500d4abSHong Zhang 3336500d4abSHong Zhang ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 3346500d4abSHong Zhang ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 3356500d4abSHong Zhang ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 3366500d4abSHong Zhang ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 3376500d4abSHong Zhang 3386500d4abSHong Zhang /* Create vertices and edges array for the subnetworks */ 3396500d4abSHong Zhang for (j=0; j < network->nsubnet; j++) { 3406500d4abSHong Zhang ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr); 3416500d4abSHong Zhang ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr); 3426500d4abSHong Zhang /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 3436500d4abSHong Zhang These get updated when the vertices and edges are added. */ 3446500d4abSHong Zhang network->subnet[j].nvtx = network->subnet[j].nedge = 0; 3456500d4abSHong Zhang } 3466500d4abSHong Zhang 3476500d4abSHong Zhang network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 3486500d4abSHong Zhang ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 3496500d4abSHong Zhang for (i=network->eStart; i < network->eEnd; i++) { 3506500d4abSHong Zhang network->header[i].index = i; /* Global edge number */ 3516500d4abSHong Zhang for (j=0; j < network->nsubnet; j++) { 3526500d4abSHong Zhang if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) { 3536500d4abSHong Zhang network->header[i].subnetid = j; /* Subnetwork id */ 3546500d4abSHong Zhang network->subnet[j].edges[network->subnet[j].nedge++] = i; 3556500d4abSHong Zhang break; 3566500d4abSHong Zhang } 3576500d4abSHong Zhang } 3586500d4abSHong Zhang network->header[i].ndata = 0; 3596500d4abSHong Zhang ndata = network->header[i].ndata; 3606500d4abSHong Zhang ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 3616500d4abSHong Zhang network->header[i].offset[ndata] = 0; 3626500d4abSHong Zhang } 3636500d4abSHong Zhang 3646500d4abSHong Zhang for(i=network->vStart; i < network->vEnd; i++) { 3656500d4abSHong Zhang network->header[i].index = i - network->vStart; 3666500d4abSHong Zhang for (j=0; j < network->nsubnet; j++) { 3676500d4abSHong Zhang if ((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) { 3686500d4abSHong Zhang network->header[i].subnetid = j; 3696500d4abSHong Zhang network->subnet[j].vertices[network->subnet[j].nvtx++] = i; 3706500d4abSHong Zhang break; 3716500d4abSHong Zhang } 3726500d4abSHong Zhang } 3736500d4abSHong Zhang network->header[i].ndata = 0; 3746500d4abSHong Zhang ndata = network->header[i].ndata; 3756500d4abSHong Zhang ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 3766500d4abSHong Zhang network->header[i].offset[ndata] = 0; 3776500d4abSHong Zhang } 3786500d4abSHong Zhang 3796500d4abSHong Zhang ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 3805f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3815f2c45f1SShri Abhyankar } 3825f2c45f1SShri Abhyankar 38394ef8ddeSSatish Balay /*@C 3842727e31bSShri Abhyankar DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork 3852727e31bSShri Abhyankar 3862727e31bSShri Abhyankar Input Parameters 3872727e31bSShri Abhyankar + dm - the number object 3882727e31bSShri Abhyankar - id - the ID (integer) of the subnetwork 3892727e31bSShri Abhyankar 3902727e31bSShri Abhyankar Output Parameters 3912727e31bSShri Abhyankar + nv - number of vertices (local) 3922727e31bSShri Abhyankar . ne - number of edges (local) 3932727e31bSShri Abhyankar . vtx - local vertices for this subnetwork 3942727e31bSShri Abhyankar . edge - local edges for this subnetwork 3952727e31bSShri Abhyankar 3962727e31bSShri Abhyankar Notes: 3972727e31bSShri Abhyankar Cannot call this routine before DMNetworkLayoutSetup() 3982727e31bSShri Abhyankar 3992727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 4002727e31bSShri Abhyankar @*/ 4012727e31bSShri Abhyankar PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge) 4022727e31bSShri Abhyankar { 4032727e31bSShri Abhyankar DM_Network *network = (DM_Network*) netdm->data; 4042727e31bSShri Abhyankar 4052727e31bSShri Abhyankar PetscFunctionBegin; 4062727e31bSShri Abhyankar *nv = network->subnet[id].nvtx; 4072727e31bSShri Abhyankar *ne = network->subnet[id].nedge; 4082727e31bSShri Abhyankar *vtx = network->subnet[id].vertices; 4092727e31bSShri Abhyankar *edge = network->subnet[id].edges; 4102727e31bSShri Abhyankar PetscFunctionReturn(0); 4112727e31bSShri Abhyankar } 4122727e31bSShri Abhyankar 4132727e31bSShri Abhyankar /*@C 4145f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 4155f2c45f1SShri Abhyankar 4165f2c45f1SShri Abhyankar Logically collective on DM 4175f2c45f1SShri Abhyankar 4185f2c45f1SShri Abhyankar Input Parameters 4195f2c45f1SShri Abhyankar + dm - the network object 4205f2c45f1SShri Abhyankar . name - the component name 4215f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 4225f2c45f1SShri Abhyankar 4235f2c45f1SShri Abhyankar Output Parameters 4245f2c45f1SShri Abhyankar . key - an integer key that defines the component 4255f2c45f1SShri Abhyankar 4265f2c45f1SShri Abhyankar Notes 4275f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 4285f2c45f1SShri Abhyankar 4295f2c45f1SShri Abhyankar Level: intermediate 4305f2c45f1SShri Abhyankar 4315f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 4325f2c45f1SShri Abhyankar @*/ 4335f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 4345f2c45f1SShri Abhyankar { 4355f2c45f1SShri Abhyankar PetscErrorCode ierr; 4365f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 4375f2c45f1SShri Abhyankar DMNetworkComponent *component=&network->component[network->ncomponent]; 4385f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 4395f2c45f1SShri Abhyankar PetscInt i; 4405f2c45f1SShri Abhyankar 4415f2c45f1SShri Abhyankar PetscFunctionBegin; 4425f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 4435f2c45f1SShri Abhyankar ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 4445f2c45f1SShri Abhyankar if (flg) { 4455f2c45f1SShri Abhyankar *key = i; 4465f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4475f2c45f1SShri Abhyankar } 4486d64e262SShri Abhyankar } 4496d64e262SShri Abhyankar if(network->ncomponent == MAX_COMPONENTS) { 4506d64e262SShri Abhyankar SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS); 4515f2c45f1SShri Abhyankar } 4525f2c45f1SShri Abhyankar 4535f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 4545f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 4555f2c45f1SShri Abhyankar *key = network->ncomponent; 4565f2c45f1SShri Abhyankar network->ncomponent++; 4575f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4585f2c45f1SShri Abhyankar } 4595f2c45f1SShri Abhyankar 4605f2c45f1SShri Abhyankar /*@ 4615f2c45f1SShri Abhyankar DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 4625f2c45f1SShri Abhyankar 4635f2c45f1SShri Abhyankar Not Collective 4645f2c45f1SShri Abhyankar 4655f2c45f1SShri Abhyankar Input Parameters: 4665f2c45f1SShri Abhyankar + dm - The DMNetwork object 4675f2c45f1SShri Abhyankar 4685f2c45f1SShri Abhyankar Output Paramters: 4695f2c45f1SShri Abhyankar + vStart - The first vertex point 4705f2c45f1SShri Abhyankar - vEnd - One beyond the last vertex point 4715f2c45f1SShri Abhyankar 4725f2c45f1SShri Abhyankar Level: intermediate 4735f2c45f1SShri Abhyankar 4745f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange 4755f2c45f1SShri Abhyankar @*/ 4765f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 4775f2c45f1SShri Abhyankar { 4785f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4795f2c45f1SShri Abhyankar 4805f2c45f1SShri Abhyankar PetscFunctionBegin; 4815f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 4825f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 4835f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4845f2c45f1SShri Abhyankar } 4855f2c45f1SShri Abhyankar 4865f2c45f1SShri Abhyankar /*@ 4875f2c45f1SShri Abhyankar DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 4885f2c45f1SShri Abhyankar 4895f2c45f1SShri Abhyankar Not Collective 4905f2c45f1SShri Abhyankar 4915f2c45f1SShri Abhyankar Input Parameters: 4925f2c45f1SShri Abhyankar + dm - The DMNetwork object 4935f2c45f1SShri Abhyankar 4945f2c45f1SShri Abhyankar Output Paramters: 4955f2c45f1SShri Abhyankar + eStart - The first edge point 4965f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 4975f2c45f1SShri Abhyankar 4985f2c45f1SShri Abhyankar Level: intermediate 4995f2c45f1SShri Abhyankar 5005f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange 5015f2c45f1SShri Abhyankar @*/ 5025f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 5035f2c45f1SShri Abhyankar { 5045f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5055f2c45f1SShri Abhyankar 5065f2c45f1SShri Abhyankar PetscFunctionBegin; 5075f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 5085f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 5095f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5105f2c45f1SShri Abhyankar } 5115f2c45f1SShri Abhyankar 5127b6afd5bSHong Zhang /*@ 513e85e6aecSHong Zhang DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge. 5147b6afd5bSHong Zhang 5157b6afd5bSHong Zhang Not Collective 5167b6afd5bSHong Zhang 5177b6afd5bSHong Zhang Input Parameters: 5187b6afd5bSHong Zhang + dm - DMNetwork object 519e85e6aecSHong Zhang - p - edge point 5207b6afd5bSHong Zhang 5217b6afd5bSHong Zhang Output Paramters: 522e85e6aecSHong Zhang . index - user global numbering for the edge 5237b6afd5bSHong Zhang 5247b6afd5bSHong Zhang Level: intermediate 5257b6afd5bSHong Zhang 526e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex 5277b6afd5bSHong Zhang @*/ 528e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index) 5297b6afd5bSHong Zhang { 5307b6afd5bSHong Zhang PetscErrorCode ierr; 5317b6afd5bSHong Zhang DM_Network *network = (DM_Network*)dm->data; 5327b6afd5bSHong Zhang PetscInt offsetp; 5337b6afd5bSHong Zhang DMNetworkComponentHeader header; 5347b6afd5bSHong Zhang 5357b6afd5bSHong Zhang PetscFunctionBegin; 5367b6afd5bSHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 5377b6afd5bSHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 538e85e6aecSHong Zhang *index = header->index; 5397b6afd5bSHong Zhang PetscFunctionReturn(0); 5407b6afd5bSHong Zhang } 5417b6afd5bSHong Zhang 5425f2c45f1SShri Abhyankar /*@ 543e85e6aecSHong Zhang DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex. 544e85e6aecSHong Zhang 545e85e6aecSHong Zhang Not Collective 546e85e6aecSHong Zhang 547e85e6aecSHong Zhang Input Parameters: 548e85e6aecSHong Zhang + dm - DMNetwork object 549e85e6aecSHong Zhang - p - vertex point 550e85e6aecSHong Zhang 551e85e6aecSHong Zhang Output Paramters: 552e85e6aecSHong Zhang . index - user global numbering for the vertex 553e85e6aecSHong Zhang 554e85e6aecSHong Zhang Level: intermediate 555e85e6aecSHong Zhang 556e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex 557e85e6aecSHong Zhang @*/ 558e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index) 559e85e6aecSHong Zhang { 560e85e6aecSHong Zhang PetscErrorCode ierr; 561e85e6aecSHong Zhang DM_Network *network = (DM_Network*)dm->data; 562e85e6aecSHong Zhang PetscInt offsetp; 563e85e6aecSHong Zhang DMNetworkComponentHeader header; 564e85e6aecSHong Zhang 565e85e6aecSHong Zhang PetscFunctionBegin; 566e85e6aecSHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 567e85e6aecSHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 568e85e6aecSHong Zhang *index = header->index; 569e85e6aecSHong Zhang PetscFunctionReturn(0); 570e85e6aecSHong Zhang } 571e85e6aecSHong Zhang 572c3b11c7cSShri Abhyankar /* 573c3b11c7cSShri Abhyankar DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the 574c3b11c7cSShri Abhyankar component value from the component data array 575c3b11c7cSShri Abhyankar 576c3b11c7cSShri Abhyankar Not Collective 577c3b11c7cSShri Abhyankar 578c3b11c7cSShri Abhyankar Input Parameters: 579c3b11c7cSShri Abhyankar + dm - The DMNetwork object 580c3b11c7cSShri Abhyankar . p - vertex/edge point 581c3b11c7cSShri Abhyankar - compnum - component number 582c3b11c7cSShri Abhyankar 583c3b11c7cSShri Abhyankar Output Parameters: 584c3b11c7cSShri Abhyankar + compkey - the key obtained when registering the component 585c3b11c7cSShri Abhyankar - offset - offset into the component data array associated with the vertex/edge point 586c3b11c7cSShri Abhyankar 587c3b11c7cSShri Abhyankar Notes: 588c3b11c7cSShri Abhyankar Typical usage: 589c3b11c7cSShri Abhyankar 590c3b11c7cSShri Abhyankar DMNetworkGetComponentDataArray(dm, &arr); 591c3b11c7cSShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 592c3b11c7cSShri Abhyankar Loop over vertices or edges 593c3b11c7cSShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 594c3b11c7cSShri Abhyankar Loop over numcomps 595c3b11c7cSShri Abhyankar DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset); 596c3b11c7cSShri Abhyankar compdata = (UserCompDataType)(arr+offset); 597c3b11c7cSShri Abhyankar 598c3b11c7cSShri Abhyankar Level: intermediate 599c3b11c7cSShri Abhyankar 600c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 601c3b11c7cSShri Abhyankar */ 602c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 603c3b11c7cSShri Abhyankar { 604c3b11c7cSShri Abhyankar PetscErrorCode ierr; 605c3b11c7cSShri Abhyankar PetscInt offsetp; 606c3b11c7cSShri Abhyankar DMNetworkComponentHeader header; 607c3b11c7cSShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 608c3b11c7cSShri Abhyankar 609c3b11c7cSShri Abhyankar PetscFunctionBegin; 610c3b11c7cSShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 611c3b11c7cSShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 612c3b11c7cSShri Abhyankar if (compkey) *compkey = header->key[compnum]; 613c3b11c7cSShri Abhyankar if (offset) *offset = offsetp+network->dataheadersize+header->offset[compnum]; 614c3b11c7cSShri Abhyankar PetscFunctionReturn(0); 615c3b11c7cSShri Abhyankar } 616c3b11c7cSShri Abhyankar 617c3b11c7cSShri Abhyankar /*@ 618c3b11c7cSShri Abhyankar DMNetworkGetComponent - Returns the network component and its key 619c3b11c7cSShri Abhyankar 620c3b11c7cSShri Abhyankar Not Collective 621c3b11c7cSShri Abhyankar 622c3b11c7cSShri Abhyankar Input Parameters 623c3b11c7cSShri Abhyankar + dm - DMNetwork object 624c3b11c7cSShri Abhyankar . p - edge or vertex point 625c3b11c7cSShri Abhyankar - compnum - component number 626c3b11c7cSShri Abhyankar 627c3b11c7cSShri Abhyankar Output Parameters: 628c3b11c7cSShri Abhyankar + compkey - the key set for this computing during registration 629c3b11c7cSShri Abhyankar - component - the component data 630c3b11c7cSShri Abhyankar 631c3b11c7cSShri Abhyankar Notes: 632c3b11c7cSShri Abhyankar Typical usage: 633c3b11c7cSShri Abhyankar 634c3b11c7cSShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 635c3b11c7cSShri Abhyankar Loop over vertices or edges 636c3b11c7cSShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 637c3b11c7cSShri Abhyankar Loop over numcomps 638c3b11c7cSShri Abhyankar DMNetworkGetComponent(dm,v,compnum,&key,&component); 639c3b11c7cSShri Abhyankar 640c3b11c7cSShri Abhyankar Level: intermediate 641c3b11c7cSShri Abhyankar 642c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset 643c3b11c7cSShri Abhyankar @*/ 644c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component) 645c3b11c7cSShri Abhyankar { 646c3b11c7cSShri Abhyankar PetscErrorCode ierr; 647c3b11c7cSShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 648c3b11c7cSShri Abhyankar PetscInt offsetd; 649c3b11c7cSShri Abhyankar 650c3b11c7cSShri Abhyankar PetscFunctionBegin; 651c3b11c7cSShri Abhyankar 652c3b11c7cSShri Abhyankar ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr); 653c3b11c7cSShri Abhyankar *component = network->componentdataarray+offsetd; 654c3b11c7cSShri Abhyankar 655c3b11c7cSShri Abhyankar PetscFunctionReturn(0); 656c3b11c7cSShri Abhyankar } 657c3b11c7cSShri Abhyankar 658e85e6aecSHong Zhang /*@ 659325661f6SSatish Balay DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 6605f2c45f1SShri Abhyankar 6615f2c45f1SShri Abhyankar Not Collective 6625f2c45f1SShri Abhyankar 6635f2c45f1SShri Abhyankar Input Parameters: 6645f2c45f1SShri Abhyankar + dm - The DMNetwork object 6655f2c45f1SShri Abhyankar . p - vertex/edge point 6665f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component 6675f2c45f1SShri Abhyankar - compvalue - pointer to the data structure for the component 6685f2c45f1SShri Abhyankar 6695f2c45f1SShri Abhyankar Level: intermediate 6705f2c45f1SShri Abhyankar 6715f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 6725f2c45f1SShri Abhyankar @*/ 6735f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 6745f2c45f1SShri Abhyankar { 6755f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 67643a39a44SBarry Smith DMNetworkComponent *component = &network->component[componentkey]; 6775f2c45f1SShri Abhyankar DMNetworkComponentHeader header = &network->header[p]; 6785f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue = &network->cvalue[p]; 6795f2c45f1SShri Abhyankar PetscErrorCode ierr; 6805f2c45f1SShri Abhyankar 6815f2c45f1SShri Abhyankar PetscFunctionBegin; 682fa58f0a9SHong 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); 683fa58f0a9SHong Zhang 68443a39a44SBarry Smith header->size[header->ndata] = component->size; 68543a39a44SBarry Smith ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 6865f2c45f1SShri Abhyankar header->key[header->ndata] = componentkey; 6875f2c45f1SShri Abhyankar if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 6885f2c45f1SShri Abhyankar 6895f2c45f1SShri Abhyankar cvalue->data[header->ndata] = (void*)compvalue; 6905f2c45f1SShri Abhyankar header->ndata++; 6915f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6925f2c45f1SShri Abhyankar } 6935f2c45f1SShri Abhyankar 6945f2c45f1SShri Abhyankar /*@ 6955f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 6965f2c45f1SShri Abhyankar 6975f2c45f1SShri Abhyankar Not Collective 6985f2c45f1SShri Abhyankar 6995f2c45f1SShri Abhyankar Input Parameters: 7005f2c45f1SShri Abhyankar + dm - The DMNetwork object 7015f2c45f1SShri Abhyankar . p - vertex/edge point 7025f2c45f1SShri Abhyankar 7035f2c45f1SShri Abhyankar Output Parameters: 7045f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 7055f2c45f1SShri Abhyankar 7065f2c45f1SShri Abhyankar Level: intermediate 7075f2c45f1SShri Abhyankar 7085f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 7095f2c45f1SShri Abhyankar @*/ 7105f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 7115f2c45f1SShri Abhyankar { 7125f2c45f1SShri Abhyankar PetscErrorCode ierr; 7135f2c45f1SShri Abhyankar PetscInt offset; 7145f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7155f2c45f1SShri Abhyankar 7165f2c45f1SShri Abhyankar PetscFunctionBegin; 7175f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 7185f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 7195f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7205f2c45f1SShri Abhyankar } 7215f2c45f1SShri Abhyankar 7225f2c45f1SShri Abhyankar /*@ 7235f2c45f1SShri Abhyankar DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 7245f2c45f1SShri Abhyankar 7255f2c45f1SShri Abhyankar Not Collective 7265f2c45f1SShri Abhyankar 7275f2c45f1SShri Abhyankar Input Parameters: 7285f2c45f1SShri Abhyankar + dm - The DMNetwork object 7295f2c45f1SShri Abhyankar - p - the edge/vertex point 7305f2c45f1SShri Abhyankar 7315f2c45f1SShri Abhyankar Output Parameters: 7325f2c45f1SShri Abhyankar . offset - the offset 7335f2c45f1SShri Abhyankar 7345f2c45f1SShri Abhyankar Level: intermediate 7355f2c45f1SShri Abhyankar 7365f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 7375f2c45f1SShri Abhyankar @*/ 7385f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 7395f2c45f1SShri Abhyankar { 7405f2c45f1SShri Abhyankar PetscErrorCode ierr; 7415f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7425f2c45f1SShri Abhyankar 7435f2c45f1SShri Abhyankar PetscFunctionBegin; 7445f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr); 7455f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7465f2c45f1SShri Abhyankar } 7475f2c45f1SShri Abhyankar 7485f2c45f1SShri Abhyankar /*@ 7495f2c45f1SShri Abhyankar DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 7505f2c45f1SShri Abhyankar 7515f2c45f1SShri Abhyankar Not Collective 7525f2c45f1SShri Abhyankar 7535f2c45f1SShri Abhyankar Input Parameters: 7545f2c45f1SShri Abhyankar + dm - The DMNetwork object 7555f2c45f1SShri Abhyankar - p - the edge/vertex point 7565f2c45f1SShri Abhyankar 7575f2c45f1SShri Abhyankar Output Parameters: 7585f2c45f1SShri Abhyankar . offsetg - the offset 7595f2c45f1SShri Abhyankar 7605f2c45f1SShri Abhyankar Level: intermediate 7615f2c45f1SShri Abhyankar 7625f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 7635f2c45f1SShri Abhyankar @*/ 7645f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 7655f2c45f1SShri Abhyankar { 7665f2c45f1SShri Abhyankar PetscErrorCode ierr; 7675f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7685f2c45f1SShri Abhyankar 7695f2c45f1SShri Abhyankar PetscFunctionBegin; 7705f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr); 7716fefedf4SHong Zhang if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */ 7725f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7735f2c45f1SShri Abhyankar } 7745f2c45f1SShri Abhyankar 77524121865SAdrian Maldonado /*@ 77624121865SAdrian Maldonado DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector. 77724121865SAdrian Maldonado 77824121865SAdrian Maldonado Not Collective 77924121865SAdrian Maldonado 78024121865SAdrian Maldonado Input Parameters: 78124121865SAdrian Maldonado + dm - The DMNetwork object 78224121865SAdrian Maldonado - p - the edge point 78324121865SAdrian Maldonado 78424121865SAdrian Maldonado Output Parameters: 78524121865SAdrian Maldonado . offset - the offset 78624121865SAdrian Maldonado 78724121865SAdrian Maldonado Level: intermediate 78824121865SAdrian Maldonado 78924121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 79024121865SAdrian Maldonado @*/ 79124121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset) 79224121865SAdrian Maldonado { 79324121865SAdrian Maldonado PetscErrorCode ierr; 79424121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 79524121865SAdrian Maldonado 79624121865SAdrian Maldonado PetscFunctionBegin; 79724121865SAdrian Maldonado 79824121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr); 79924121865SAdrian Maldonado PetscFunctionReturn(0); 80024121865SAdrian Maldonado } 80124121865SAdrian Maldonado 80224121865SAdrian Maldonado /*@ 80324121865SAdrian Maldonado DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector. 80424121865SAdrian Maldonado 80524121865SAdrian Maldonado Not Collective 80624121865SAdrian Maldonado 80724121865SAdrian Maldonado Input Parameters: 80824121865SAdrian Maldonado + dm - The DMNetwork object 80924121865SAdrian Maldonado - p - the vertex point 81024121865SAdrian Maldonado 81124121865SAdrian Maldonado Output Parameters: 81224121865SAdrian Maldonado . offset - the offset 81324121865SAdrian Maldonado 81424121865SAdrian Maldonado Level: intermediate 81524121865SAdrian Maldonado 81624121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 81724121865SAdrian Maldonado @*/ 81824121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset) 81924121865SAdrian Maldonado { 82024121865SAdrian Maldonado PetscErrorCode ierr; 82124121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 82224121865SAdrian Maldonado 82324121865SAdrian Maldonado PetscFunctionBegin; 82424121865SAdrian Maldonado 82524121865SAdrian Maldonado p -= network->vStart; 82624121865SAdrian Maldonado 82724121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr); 82824121865SAdrian Maldonado PetscFunctionReturn(0); 82924121865SAdrian Maldonado } 8305f2c45f1SShri Abhyankar /*@ 8315f2c45f1SShri Abhyankar DMNetworkAddNumVariables - Add number of variables associated with a given point. 8325f2c45f1SShri Abhyankar 8335f2c45f1SShri Abhyankar Not Collective 8345f2c45f1SShri Abhyankar 8355f2c45f1SShri Abhyankar Input Parameters: 8365f2c45f1SShri Abhyankar + dm - The DMNetworkObject 8375f2c45f1SShri Abhyankar . p - the vertex/edge point 8385f2c45f1SShri Abhyankar - nvar - number of additional variables 8395f2c45f1SShri Abhyankar 8405f2c45f1SShri Abhyankar Level: intermediate 8415f2c45f1SShri Abhyankar 8425f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables 8435f2c45f1SShri Abhyankar @*/ 8445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 8455f2c45f1SShri Abhyankar { 8465f2c45f1SShri Abhyankar PetscErrorCode ierr; 8475f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8485f2c45f1SShri Abhyankar 8495f2c45f1SShri Abhyankar PetscFunctionBegin; 8505f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 8515f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8525f2c45f1SShri Abhyankar } 8535f2c45f1SShri Abhyankar 85427f51fceSHong Zhang /*@ 85527f51fceSHong Zhang DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point. 85627f51fceSHong Zhang 85727f51fceSHong Zhang Not Collective 85827f51fceSHong Zhang 85927f51fceSHong Zhang Input Parameters: 86027f51fceSHong Zhang + dm - The DMNetworkObject 86127f51fceSHong Zhang - p - the vertex/edge point 86227f51fceSHong Zhang 86327f51fceSHong Zhang Output Parameters: 86427f51fceSHong Zhang . nvar - number of variables 86527f51fceSHong Zhang 86627f51fceSHong Zhang Level: intermediate 86727f51fceSHong Zhang 86827f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables 86927f51fceSHong Zhang @*/ 87027f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar) 87127f51fceSHong Zhang { 87227f51fceSHong Zhang PetscErrorCode ierr; 87327f51fceSHong Zhang DM_Network *network = (DM_Network*)dm->data; 87427f51fceSHong Zhang 87527f51fceSHong Zhang PetscFunctionBegin; 87627f51fceSHong Zhang ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 87727f51fceSHong Zhang PetscFunctionReturn(0); 87827f51fceSHong Zhang } 87927f51fceSHong Zhang 8805f2c45f1SShri Abhyankar /*@ 8815f2c45f1SShri Abhyankar DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 8825f2c45f1SShri Abhyankar 8835f2c45f1SShri Abhyankar Not Collective 8845f2c45f1SShri Abhyankar 8855f2c45f1SShri Abhyankar Input Parameters: 8865f2c45f1SShri Abhyankar + dm - The DMNetworkObject 8875f2c45f1SShri Abhyankar . p - the vertex/edge point 8885f2c45f1SShri Abhyankar - nvar - number of variables 8895f2c45f1SShri Abhyankar 8905f2c45f1SShri Abhyankar Level: intermediate 8915f2c45f1SShri Abhyankar 8925f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables 8935f2c45f1SShri Abhyankar @*/ 8945f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 8955f2c45f1SShri Abhyankar { 8965f2c45f1SShri Abhyankar PetscErrorCode ierr; 8975f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8985f2c45f1SShri Abhyankar 8995f2c45f1SShri Abhyankar PetscFunctionBegin; 9005f2c45f1SShri Abhyankar ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 9015f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9025f2c45f1SShri Abhyankar } 9035f2c45f1SShri Abhyankar 9045f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This 9055f2c45f1SShri Abhyankar function is called during DMSetUp() */ 9065f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 9075f2c45f1SShri Abhyankar { 9085f2c45f1SShri Abhyankar PetscErrorCode ierr; 9095f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9105f2c45f1SShri Abhyankar PetscInt arr_size; 9115f2c45f1SShri Abhyankar PetscInt p,offset,offsetp; 9125f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 9135f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 9145f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 9155f2c45f1SShri Abhyankar PetscInt ncomp, i; 9165f2c45f1SShri Abhyankar 9175f2c45f1SShri Abhyankar PetscFunctionBegin; 9185f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 9195f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 92075b160a0SShri Abhyankar ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 9215f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 9225f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 9235f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 9245f2c45f1SShri Abhyankar /* Copy header */ 9255f2c45f1SShri Abhyankar header = &network->header[p]; 926302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 9275f2c45f1SShri Abhyankar /* Copy data */ 9285f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 9295f2c45f1SShri Abhyankar ncomp = header->ndata; 9305f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 9315f2c45f1SShri Abhyankar offset = offsetp + network->dataheadersize + header->offset[i]; 932302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 9335f2c45f1SShri Abhyankar } 9345f2c45f1SShri Abhyankar } 9355f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9365f2c45f1SShri Abhyankar } 9375f2c45f1SShri Abhyankar 9385f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 9395f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm) 9405f2c45f1SShri Abhyankar { 9415f2c45f1SShri Abhyankar PetscErrorCode ierr; 9425f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9435f2c45f1SShri Abhyankar 9445f2c45f1SShri Abhyankar PetscFunctionBegin; 9455f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 9465f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9475f2c45f1SShri Abhyankar } 9485f2c45f1SShri Abhyankar 9495f2c45f1SShri Abhyankar /*@C 9505f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray - Returns the component data array 9515f2c45f1SShri Abhyankar 9525f2c45f1SShri Abhyankar Not Collective 9535f2c45f1SShri Abhyankar 9545f2c45f1SShri Abhyankar Input Parameters: 9555f2c45f1SShri Abhyankar . dm - The DMNetwork Object 9565f2c45f1SShri Abhyankar 9575f2c45f1SShri Abhyankar Output Parameters: 9585f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components 9595f2c45f1SShri Abhyankar 9605f2c45f1SShri Abhyankar Level: intermediate 9615f2c45f1SShri Abhyankar 962a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents 9635f2c45f1SShri Abhyankar @*/ 9645f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 9655f2c45f1SShri Abhyankar { 9665f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9675f2c45f1SShri Abhyankar 9685f2c45f1SShri Abhyankar PetscFunctionBegin; 9695f2c45f1SShri Abhyankar *componentdataarray = network->componentdataarray; 9705f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9715f2c45f1SShri Abhyankar } 9725f2c45f1SShri Abhyankar 97324121865SAdrian Maldonado /* Get a subsection from a range of points */ 97424121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection) 97524121865SAdrian Maldonado { 97624121865SAdrian Maldonado PetscErrorCode ierr; 97724121865SAdrian Maldonado PetscInt i, nvar; 97824121865SAdrian Maldonado 97924121865SAdrian Maldonado PetscFunctionBegin; 98024121865SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr); 98124121865SAdrian Maldonado ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr); 98224121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 98324121865SAdrian Maldonado ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr); 98424121865SAdrian Maldonado ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr); 98524121865SAdrian Maldonado } 98624121865SAdrian Maldonado 98724121865SAdrian Maldonado ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr); 98824121865SAdrian Maldonado PetscFunctionReturn(0); 98924121865SAdrian Maldonado } 99024121865SAdrian Maldonado 99124121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */ 99224121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 99324121865SAdrian Maldonado { 99424121865SAdrian Maldonado PetscErrorCode ierr; 99524121865SAdrian Maldonado PetscInt i, *subpoints; 99624121865SAdrian Maldonado 99724121865SAdrian Maldonado PetscFunctionBegin; 99824121865SAdrian Maldonado /* Create index sets to map from "points" to "subpoints" */ 99924121865SAdrian Maldonado ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr); 100024121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 100124121865SAdrian Maldonado subpoints[i - pstart] = i; 100224121865SAdrian Maldonado } 1003459726d8SSatish Balay ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr); 100424121865SAdrian Maldonado ierr = PetscFree(subpoints);CHKERRQ(ierr); 100524121865SAdrian Maldonado PetscFunctionReturn(0); 100624121865SAdrian Maldonado } 100724121865SAdrian Maldonado 100824121865SAdrian Maldonado /*@ 100924121865SAdrian Maldonado DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute. 101024121865SAdrian Maldonado 101124121865SAdrian Maldonado Collective 101224121865SAdrian Maldonado 101324121865SAdrian Maldonado Input Parameters: 101424121865SAdrian Maldonado . dm - The DMNetworkObject 101524121865SAdrian Maldonado 101624121865SAdrian Maldonado Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 101724121865SAdrian Maldonado 101824121865SAdrian Maldonado points = [0 1 2 3 4 5 6] 101924121865SAdrian Maldonado 102024121865SAdrian Maldonado where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]). 102124121865SAdrian Maldonado 102224121865SAdrian Maldonado With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 102324121865SAdrian Maldonado 102424121865SAdrian Maldonado Level: intermediate 102524121865SAdrian Maldonado 102624121865SAdrian Maldonado @*/ 102724121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 102824121865SAdrian Maldonado { 102924121865SAdrian Maldonado PetscErrorCode ierr; 103024121865SAdrian Maldonado MPI_Comm comm; 10319852e123SBarry Smith PetscMPIInt rank, size; 103224121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 103324121865SAdrian Maldonado 1034eab1376dSHong Zhang PetscFunctionBegin; 103524121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 103624121865SAdrian Maldonado ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 10379852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 103824121865SAdrian Maldonado 103924121865SAdrian Maldonado /* Create maps for vertices and edges */ 104024121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 104124121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr); 104224121865SAdrian Maldonado 104324121865SAdrian Maldonado /* Create local sub-sections */ 104424121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr); 104524121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr); 104624121865SAdrian Maldonado 10479852e123SBarry Smith if (size > 1) { 104824121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 104924121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr); 105024121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr); 105124121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr); 105224121865SAdrian Maldonado } else { 105324121865SAdrian Maldonado /* create structures for vertex */ 105424121865SAdrian Maldonado ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr); 105524121865SAdrian Maldonado /* create structures for edge */ 105624121865SAdrian Maldonado ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr); 105724121865SAdrian Maldonado } 105824121865SAdrian Maldonado 105924121865SAdrian Maldonado 106024121865SAdrian Maldonado /* Add viewers */ 106124121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr); 106224121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr); 106324121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr); 106424121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr); 106524121865SAdrian Maldonado 106624121865SAdrian Maldonado PetscFunctionReturn(0); 106724121865SAdrian Maldonado } 10687b6afd5bSHong Zhang 10695f2c45f1SShri Abhyankar /*@ 10705f2c45f1SShri Abhyankar DMNetworkDistribute - Distributes the network and moves associated component data. 10715f2c45f1SShri Abhyankar 10725f2c45f1SShri Abhyankar Collective 10735f2c45f1SShri Abhyankar 10745f2c45f1SShri Abhyankar Input Parameter: 1075d3464fd4SAdrian Maldonado + DM - the DMNetwork object 10765f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default 10775f2c45f1SShri Abhyankar 10785f2c45f1SShri Abhyankar Notes: 10798b171c8eSHong Zhang Distributes the network with <overlap>-overlapping partitioning of the edges. 10805f2c45f1SShri Abhyankar 10815f2c45f1SShri Abhyankar Level: intermediate 10825f2c45f1SShri Abhyankar 10835f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 10845f2c45f1SShri Abhyankar @*/ 1085d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap) 10865f2c45f1SShri Abhyankar { 1087d3464fd4SAdrian Maldonado MPI_Comm comm; 10885f2c45f1SShri Abhyankar PetscErrorCode ierr; 1089d3464fd4SAdrian Maldonado PetscMPIInt size; 1090d3464fd4SAdrian Maldonado DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data); 1091d3464fd4SAdrian Maldonado DM_Network *newDMnetwork; 10925f2c45f1SShri Abhyankar PetscSF pointsf; 10935f2c45f1SShri Abhyankar DM newDM; 109451ac5effSHong Zhang PetscPartitioner part; 1095b9c6e19dSShri Abhyankar PetscInt j,e,v,offset; 1096b9c6e19dSShri Abhyankar DMNetworkComponentHeader header; 10975f2c45f1SShri Abhyankar 10985f2c45f1SShri Abhyankar PetscFunctionBegin; 1099d3464fd4SAdrian Maldonado 1100d3464fd4SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr); 1101d3464fd4SAdrian Maldonado ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1102d3464fd4SAdrian Maldonado if (size == 1) PetscFunctionReturn(0); 1103d3464fd4SAdrian Maldonado 1104d3464fd4SAdrian Maldonado ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr); 11055f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 11065f2c45f1SShri Abhyankar newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 110751ac5effSHong Zhang 110851ac5effSHong Zhang /* Enable runtime options for petscpartitioner */ 110951ac5effSHong Zhang ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr); 111051ac5effSHong Zhang ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 111151ac5effSHong Zhang 11125f2c45f1SShri Abhyankar /* Distribute plex dm and dof section */ 111380cf41d5SMatthew G. Knepley ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 111451ac5effSHong Zhang 11155f2c45f1SShri Abhyankar /* Distribute dof section */ 1116d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr); 11175f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 1118d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr); 111951ac5effSHong Zhang 11205f2c45f1SShri Abhyankar /* Distribute data and associated section */ 112131da1fc8SHong Zhang ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 112224121865SAdrian Maldonado 11235f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 11245f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 11255f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 11265f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 11276fefedf4SHong Zhang newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart; 11286fefedf4SHong Zhang newDMnetwork->NVertices = oldDMnetwork->NVertices; 11295f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 113024121865SAdrian Maldonado 11315f2c45f1SShri Abhyankar /* Set Dof section as the default section for dm */ 11325f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 11335f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 11345f2c45f1SShri Abhyankar 1135b9c6e19dSShri Abhyankar /* Set up subnetwork info in the newDM */ 1136b9c6e19dSShri Abhyankar newDMnetwork->nsubnet = oldDMnetwork->nsubnet; 1137b9c6e19dSShri Abhyankar ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr); 1138b9c6e19dSShri Abhyankar /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already 1139b9c6e19dSShri Abhyankar calculated in DMNetworkLayoutSetUp() 1140b9c6e19dSShri Abhyankar */ 1141b9c6e19dSShri Abhyankar for(j=0; j < newDMnetwork->nsubnet; j++) { 1142b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx; 1143b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge; 1144b9c6e19dSShri Abhyankar } 1145b9c6e19dSShri Abhyankar 1146b9c6e19dSShri Abhyankar for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) { 1147b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1148b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1149b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].nedge++; 1150b9c6e19dSShri Abhyankar } 1151b9c6e19dSShri Abhyankar 1152b9c6e19dSShri Abhyankar for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) { 1153b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1154b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1155b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].nvtx++; 1156b9c6e19dSShri Abhyankar } 1157b9c6e19dSShri Abhyankar 1158b9c6e19dSShri Abhyankar /* Now create the vertices and edge arrays for the subnetworks */ 1159b9c6e19dSShri Abhyankar for(j=0; j < newDMnetwork->nsubnet; j++) { 1160b9c6e19dSShri Abhyankar ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr); 1161b9c6e19dSShri Abhyankar ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr); 1162b9c6e19dSShri Abhyankar /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 1163b9c6e19dSShri Abhyankar These get updated when the vertices and edges are added. */ 1164b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0; 1165b9c6e19dSShri Abhyankar } 1166b9c6e19dSShri Abhyankar 1167b9c6e19dSShri Abhyankar /* Set the vertices and edges in each subnetwork */ 1168b9c6e19dSShri Abhyankar for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) { 1169b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1170b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1171b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e; 1172b9c6e19dSShri Abhyankar } 1173b9c6e19dSShri Abhyankar 1174b9c6e19dSShri Abhyankar for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) { 1175b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1176b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1177b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v; 1178b9c6e19dSShri Abhyankar } 1179b9c6e19dSShri Abhyankar 118024121865SAdrian Maldonado /* Destroy point SF */ 118124121865SAdrian Maldonado ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 118224121865SAdrian Maldonado 1183d3464fd4SAdrian Maldonado ierr = DMDestroy(dm);CHKERRQ(ierr); 1184d3464fd4SAdrian Maldonado *dm = newDM; 11855f2c45f1SShri Abhyankar PetscFunctionReturn(0); 11865f2c45f1SShri Abhyankar } 11875f2c45f1SShri Abhyankar 118824121865SAdrian Maldonado /*@C 118924121865SAdrian Maldonado PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering. 119024121865SAdrian Maldonado 119124121865SAdrian Maldonado Input Parameters: 119224121865SAdrian Maldonado + masterSF - the original SF structure 119324121865SAdrian Maldonado - map - a ISLocalToGlobal mapping that contains the subset of points 119424121865SAdrian Maldonado 119524121865SAdrian Maldonado Output Parameters: 119624121865SAdrian Maldonado . subSF - a subset of the masterSF for the desired subset. 119724121865SAdrian Maldonado */ 119824121865SAdrian Maldonado 119924121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) { 120024121865SAdrian Maldonado 120124121865SAdrian Maldonado PetscErrorCode ierr; 120224121865SAdrian Maldonado PetscInt nroots, nleaves, *ilocal_sub; 120324121865SAdrian Maldonado PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 120424121865SAdrian Maldonado PetscInt *local_points, *remote_points; 120524121865SAdrian Maldonado PetscSFNode *iremote_sub; 120624121865SAdrian Maldonado const PetscInt *ilocal; 120724121865SAdrian Maldonado const PetscSFNode *iremote; 120824121865SAdrian Maldonado 120924121865SAdrian Maldonado PetscFunctionBegin; 121024121865SAdrian Maldonado ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 121124121865SAdrian Maldonado 121224121865SAdrian Maldonado /* Look for leaves that pertain to the subset of points. Get the local ordering */ 121324121865SAdrian Maldonado ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr); 121424121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr); 121524121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 121624121865SAdrian Maldonado if (ilocal_map[i] != -1) nleaves_sub += 1; 121724121865SAdrian Maldonado } 121824121865SAdrian Maldonado /* Re-number ilocal with subset numbering. Need information from roots */ 121924121865SAdrian Maldonado ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr); 122024121865SAdrian Maldonado for (i = 0; i < nroots; i++) local_points[i] = i; 122124121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr); 122224121865SAdrian Maldonado ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 122324121865SAdrian Maldonado ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 122424121865SAdrian Maldonado /* Fill up graph using local (that is, local to the subset) numbering. */ 12254b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr); 12264b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr); 122724121865SAdrian Maldonado nleaves_sub = 0; 122824121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 122924121865SAdrian Maldonado if (ilocal_map[i] != -1) { 123024121865SAdrian Maldonado ilocal_sub[nleaves_sub] = ilocal_map[i]; 12314b70a8deSAdrian Maldonado iremote_sub[nleaves_sub].rank = iremote[i].rank; 123224121865SAdrian Maldonado iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 123324121865SAdrian Maldonado nleaves_sub += 1; 123424121865SAdrian Maldonado } 123524121865SAdrian Maldonado } 123624121865SAdrian Maldonado ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr); 123724121865SAdrian Maldonado ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr); 123824121865SAdrian Maldonado 123924121865SAdrian Maldonado /* Create new subSF */ 124024121865SAdrian Maldonado ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr); 124124121865SAdrian Maldonado ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr); 12424b70a8deSAdrian Maldonado ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr); 124324121865SAdrian Maldonado ierr = PetscFree(ilocal_map);CHKERRQ(ierr); 12444b70a8deSAdrian Maldonado ierr = PetscFree(iremote_sub);CHKERRQ(ierr); 124524121865SAdrian Maldonado PetscFunctionReturn(0); 124624121865SAdrian Maldonado } 124724121865SAdrian Maldonado 12485f2c45f1SShri Abhyankar /*@C 12495f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 12505f2c45f1SShri Abhyankar 12515f2c45f1SShri Abhyankar Not Collective 12525f2c45f1SShri Abhyankar 12535f2c45f1SShri Abhyankar Input Parameters: 12545f2c45f1SShri Abhyankar + dm - The DMNetwork object 12555f2c45f1SShri Abhyankar - p - the vertex point 12565f2c45f1SShri Abhyankar 12575f2c45f1SShri Abhyankar Output Paramters: 12585f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 12595f2c45f1SShri Abhyankar - edges - List of edge points 12605f2c45f1SShri Abhyankar 12615f2c45f1SShri Abhyankar Level: intermediate 12625f2c45f1SShri Abhyankar 12635f2c45f1SShri Abhyankar Fortran Notes: 12645f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 12655f2c45f1SShri Abhyankar include petsc.h90 in your code. 12665f2c45f1SShri Abhyankar 1267d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices 12685f2c45f1SShri Abhyankar @*/ 12695f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 12705f2c45f1SShri Abhyankar { 12715f2c45f1SShri Abhyankar PetscErrorCode ierr; 12725f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 12735f2c45f1SShri Abhyankar 12745f2c45f1SShri Abhyankar PetscFunctionBegin; 12755f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 12765f2c45f1SShri Abhyankar ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 12775f2c45f1SShri Abhyankar PetscFunctionReturn(0); 12785f2c45f1SShri Abhyankar } 12795f2c45f1SShri Abhyankar 12805f2c45f1SShri Abhyankar /*@C 1281d842c372SHong Zhang DMNetworkGetConnectedVertices - Return the connected vertices for this edge point 12825f2c45f1SShri Abhyankar 12835f2c45f1SShri Abhyankar Not Collective 12845f2c45f1SShri Abhyankar 12855f2c45f1SShri Abhyankar Input Parameters: 12865f2c45f1SShri Abhyankar + dm - The DMNetwork object 12875f2c45f1SShri Abhyankar - p - the edge point 12885f2c45f1SShri Abhyankar 12895f2c45f1SShri Abhyankar Output Paramters: 12905f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 12915f2c45f1SShri Abhyankar 12925f2c45f1SShri Abhyankar Level: intermediate 12935f2c45f1SShri Abhyankar 12945f2c45f1SShri Abhyankar Fortran Notes: 12955f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 12965f2c45f1SShri Abhyankar include petsc.h90 in your code. 12975f2c45f1SShri Abhyankar 12985f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 12995f2c45f1SShri Abhyankar @*/ 1300d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[]) 13015f2c45f1SShri Abhyankar { 13025f2c45f1SShri Abhyankar PetscErrorCode ierr; 13035f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 13045f2c45f1SShri Abhyankar 13055f2c45f1SShri Abhyankar PetscFunctionBegin; 13065f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 13075f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13085f2c45f1SShri Abhyankar } 13095f2c45f1SShri Abhyankar 13105f2c45f1SShri Abhyankar /*@ 13115f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 13125f2c45f1SShri Abhyankar 13135f2c45f1SShri Abhyankar Not Collective 13145f2c45f1SShri Abhyankar 13155f2c45f1SShri Abhyankar Input Parameters: 13165f2c45f1SShri Abhyankar + dm - The DMNetwork object 13175f2c45f1SShri Abhyankar . p - the vertex point 13185f2c45f1SShri Abhyankar 13195f2c45f1SShri Abhyankar Output Parameter: 13205f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 13215f2c45f1SShri Abhyankar 13225f2c45f1SShri Abhyankar Level: intermediate 13235f2c45f1SShri Abhyankar 1324d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange 13255f2c45f1SShri Abhyankar @*/ 13265f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 13275f2c45f1SShri Abhyankar { 13285f2c45f1SShri Abhyankar PetscErrorCode ierr; 13295f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 13305f2c45f1SShri Abhyankar PetscInt offsetg; 13315f2c45f1SShri Abhyankar PetscSection sectiong; 13325f2c45f1SShri Abhyankar 13335f2c45f1SShri Abhyankar PetscFunctionBegin; 13345f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 13355f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 13365f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 13375f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 13385f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13395f2c45f1SShri Abhyankar } 13405f2c45f1SShri Abhyankar 13415f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 13425f2c45f1SShri Abhyankar { 13435f2c45f1SShri Abhyankar PetscErrorCode ierr; 13445f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 13455f2c45f1SShri Abhyankar 13465f2c45f1SShri Abhyankar PetscFunctionBegin; 13475f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 13485f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 13495f2c45f1SShri Abhyankar 13505f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 13515f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 13525f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13535f2c45f1SShri Abhyankar } 13545f2c45f1SShri Abhyankar 13551ad426b7SHong Zhang /*@ 135617df6e9eSHong Zhang DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 13571ad426b7SHong Zhang -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 13581ad426b7SHong Zhang 13591ad426b7SHong Zhang Collective 13601ad426b7SHong Zhang 13611ad426b7SHong Zhang Input Parameters: 136283b2e829SHong Zhang + dm - The DMNetwork object 136383b2e829SHong Zhang . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 136483b2e829SHong Zhang - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 13651ad426b7SHong Zhang 13661ad426b7SHong Zhang Level: intermediate 13671ad426b7SHong Zhang 13681ad426b7SHong Zhang @*/ 136983b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 13701ad426b7SHong Zhang { 13711ad426b7SHong Zhang DM_Network *network=(DM_Network*)dm->data; 13728675203cSHong Zhang PetscErrorCode ierr; 13731ad426b7SHong Zhang 13741ad426b7SHong Zhang PetscFunctionBegin; 137583b2e829SHong Zhang network->userEdgeJacobian = eflg; 137683b2e829SHong Zhang network->userVertexJacobian = vflg; 13778675203cSHong Zhang 13788675203cSHong Zhang if (eflg && !network->Je) { 13798675203cSHong Zhang ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 13808675203cSHong Zhang } 13818675203cSHong Zhang 13828675203cSHong Zhang if (vflg && !network->Jv) { 13838675203cSHong Zhang PetscInt i,*vptr,nedges,vStart=network->vStart; 13848675203cSHong Zhang PetscInt nVertices = network->nVertices,nedges_total; 13858675203cSHong Zhang const PetscInt *edges; 13868675203cSHong Zhang 13878675203cSHong Zhang /* count nvertex_total */ 13888675203cSHong Zhang nedges_total = 0; 13898675203cSHong Zhang ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr); 13908675203cSHong Zhang 13918675203cSHong Zhang vptr[0] = 0; 13928675203cSHong Zhang for (i=0; i<nVertices; i++) { 13938675203cSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 13948675203cSHong Zhang nedges_total += nedges; 13958675203cSHong Zhang vptr[i+1] = vptr[i] + 2*nedges + 1; 13968675203cSHong Zhang } 13978675203cSHong Zhang 13988675203cSHong Zhang ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr); 13998675203cSHong Zhang network->Jvptr = vptr; 14008675203cSHong Zhang } 14011ad426b7SHong Zhang PetscFunctionReturn(0); 14021ad426b7SHong Zhang } 14031ad426b7SHong Zhang 14041ad426b7SHong Zhang /*@ 140583b2e829SHong Zhang DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 140683b2e829SHong Zhang 140783b2e829SHong Zhang Not Collective 140883b2e829SHong Zhang 140983b2e829SHong Zhang Input Parameters: 141083b2e829SHong Zhang + dm - The DMNetwork object 141183b2e829SHong Zhang . p - the edge point 14123e97b6e8SHong Zhang - J - array (size = 3) of Jacobian submatrices for this edge point: 14133e97b6e8SHong Zhang J[0]: this edge 1414d842c372SHong Zhang J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices() 141583b2e829SHong Zhang 141683b2e829SHong Zhang Level: intermediate 141783b2e829SHong Zhang 141883b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix 141983b2e829SHong Zhang @*/ 142083b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 142183b2e829SHong Zhang { 142283b2e829SHong Zhang DM_Network *network=(DM_Network*)dm->data; 142383b2e829SHong Zhang 142483b2e829SHong Zhang PetscFunctionBegin; 14258675203cSHong Zhang if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 14268675203cSHong Zhang 14278675203cSHong Zhang if (J) { 1428883e35e8SHong Zhang network->Je[3*p] = J[0]; 1429883e35e8SHong Zhang network->Je[3*p+1] = J[1]; 1430883e35e8SHong Zhang network->Je[3*p+2] = J[2]; 14318675203cSHong Zhang } 143283b2e829SHong Zhang PetscFunctionReturn(0); 143383b2e829SHong Zhang } 143483b2e829SHong Zhang 143583b2e829SHong Zhang /*@ 143676ddfea5SHong Zhang DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 14371ad426b7SHong Zhang 14381ad426b7SHong Zhang Not Collective 14391ad426b7SHong Zhang 14401ad426b7SHong Zhang Input Parameters: 14411ad426b7SHong Zhang + dm - The DMNetwork object 14421ad426b7SHong Zhang . p - the vertex point 14433e97b6e8SHong Zhang - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 14443e97b6e8SHong Zhang J[0]: this vertex 14453e97b6e8SHong Zhang J[1+2*i]: i-th supporting edge 14463e97b6e8SHong Zhang J[1+2*i+1]: i-th connected vertex 14471ad426b7SHong Zhang 14481ad426b7SHong Zhang Level: intermediate 14491ad426b7SHong Zhang 145083b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix 14511ad426b7SHong Zhang @*/ 1452883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 14535f2c45f1SShri Abhyankar { 14545f2c45f1SShri Abhyankar PetscErrorCode ierr; 14555f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 14568675203cSHong Zhang PetscInt i,*vptr,nedges,vStart=network->vStart; 1457883e35e8SHong Zhang const PetscInt *edges; 14585f2c45f1SShri Abhyankar 14595f2c45f1SShri Abhyankar PetscFunctionBegin; 14608675203cSHong Zhang if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 1461883e35e8SHong Zhang 14628675203cSHong Zhang if (J) { 1463883e35e8SHong Zhang vptr = network->Jvptr; 14643e97b6e8SHong Zhang network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 14653e97b6e8SHong Zhang 14663e97b6e8SHong Zhang /* Set Jacobian for each supporting edge and connected vertex */ 1467883e35e8SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 1468883e35e8SHong Zhang for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 14698675203cSHong Zhang } 1470883e35e8SHong Zhang PetscFunctionReturn(0); 1471883e35e8SHong Zhang } 1472883e35e8SHong Zhang 1473e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 14745cf7da58SHong Zhang { 14755cf7da58SHong Zhang PetscErrorCode ierr; 14765cf7da58SHong Zhang PetscInt j; 14775cf7da58SHong Zhang PetscScalar val=(PetscScalar)ncols; 14785cf7da58SHong Zhang 14795cf7da58SHong Zhang PetscFunctionBegin; 14805cf7da58SHong Zhang if (!ghost) { 14815cf7da58SHong Zhang for (j=0; j<nrows; j++) { 14825cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 14835cf7da58SHong Zhang } 14845cf7da58SHong Zhang } else { 14855cf7da58SHong Zhang for (j=0; j<nrows; j++) { 14865cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 14875cf7da58SHong Zhang } 14885cf7da58SHong Zhang } 14895cf7da58SHong Zhang PetscFunctionReturn(0); 14905cf7da58SHong Zhang } 14915cf7da58SHong Zhang 1492e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 14935cf7da58SHong Zhang { 14945cf7da58SHong Zhang PetscErrorCode ierr; 14955cf7da58SHong Zhang PetscInt j,ncols_u; 14965cf7da58SHong Zhang PetscScalar val; 14975cf7da58SHong Zhang 14985cf7da58SHong Zhang PetscFunctionBegin; 14995cf7da58SHong Zhang if (!ghost) { 15005cf7da58SHong Zhang for (j=0; j<nrows; j++) { 15015cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 15025cf7da58SHong Zhang val = (PetscScalar)ncols_u; 15035cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 15045cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 15055cf7da58SHong Zhang } 15065cf7da58SHong Zhang } else { 15075cf7da58SHong Zhang for (j=0; j<nrows; j++) { 15085cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 15095cf7da58SHong Zhang val = (PetscScalar)ncols_u; 15105cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 15115cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 15125cf7da58SHong Zhang } 15135cf7da58SHong Zhang } 15145cf7da58SHong Zhang PetscFunctionReturn(0); 15155cf7da58SHong Zhang } 15165cf7da58SHong Zhang 1517e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 15185cf7da58SHong Zhang { 15195cf7da58SHong Zhang PetscErrorCode ierr; 15205cf7da58SHong Zhang 15215cf7da58SHong Zhang PetscFunctionBegin; 15225cf7da58SHong Zhang if (Ju) { 15235cf7da58SHong Zhang ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 15245cf7da58SHong Zhang } else { 15255cf7da58SHong Zhang ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 15265cf7da58SHong Zhang } 15275cf7da58SHong Zhang PetscFunctionReturn(0); 15285cf7da58SHong Zhang } 15295cf7da58SHong Zhang 1530e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1531883e35e8SHong Zhang { 1532883e35e8SHong Zhang PetscErrorCode ierr; 1533883e35e8SHong Zhang PetscInt j,*cols; 1534883e35e8SHong Zhang PetscScalar *zeros; 1535883e35e8SHong Zhang 1536883e35e8SHong Zhang PetscFunctionBegin; 1537883e35e8SHong Zhang ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 1538883e35e8SHong Zhang for (j=0; j<ncols; j++) cols[j] = j+ cstart; 1539883e35e8SHong Zhang ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 1540883e35e8SHong Zhang ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 15411ad426b7SHong Zhang PetscFunctionReturn(0); 15421ad426b7SHong Zhang } 1543a4e85ca8SHong Zhang 1544e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 15453e97b6e8SHong Zhang { 15463e97b6e8SHong Zhang PetscErrorCode ierr; 15473e97b6e8SHong Zhang PetscInt j,M,N,row,col,ncols_u; 15483e97b6e8SHong Zhang const PetscInt *cols; 15493e97b6e8SHong Zhang PetscScalar zero=0.0; 15503e97b6e8SHong Zhang 15513e97b6e8SHong Zhang PetscFunctionBegin; 15523e97b6e8SHong Zhang ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 15533e97b6e8SHong 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); 15543e97b6e8SHong Zhang 15553e97b6e8SHong Zhang for (row=0; row<nrows; row++) { 15563e97b6e8SHong Zhang ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 15573e97b6e8SHong Zhang for (j=0; j<ncols_u; j++) { 15583e97b6e8SHong Zhang col = cols[j] + cstart; 15593e97b6e8SHong Zhang ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 15603e97b6e8SHong Zhang } 15613e97b6e8SHong Zhang ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 15623e97b6e8SHong Zhang } 15633e97b6e8SHong Zhang PetscFunctionReturn(0); 15643e97b6e8SHong Zhang } 15651ad426b7SHong Zhang 1566e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1567a4e85ca8SHong Zhang { 1568a4e85ca8SHong Zhang PetscErrorCode ierr; 1569f4431b8cSHong Zhang 1570a4e85ca8SHong Zhang PetscFunctionBegin; 1571a4e85ca8SHong Zhang if (Ju) { 1572a4e85ca8SHong Zhang ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1573a4e85ca8SHong Zhang } else { 1574a4e85ca8SHong Zhang ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1575a4e85ca8SHong Zhang } 1576a4e85ca8SHong Zhang PetscFunctionReturn(0); 1577a4e85ca8SHong Zhang } 1578a4e85ca8SHong Zhang 157924121865SAdrian 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. 158024121865SAdrian Maldonado */ 158124121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 158224121865SAdrian Maldonado { 158324121865SAdrian Maldonado PetscErrorCode ierr; 158424121865SAdrian Maldonado PetscInt i, size, dof; 158524121865SAdrian Maldonado PetscInt *glob2loc; 158624121865SAdrian Maldonado 158724121865SAdrian Maldonado PetscFunctionBegin; 158824121865SAdrian Maldonado ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr); 158924121865SAdrian Maldonado ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr); 159024121865SAdrian Maldonado 159124121865SAdrian Maldonado for (i = 0; i < size; i++) { 159224121865SAdrian Maldonado ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr); 159324121865SAdrian Maldonado dof = (dof >= 0) ? dof : -(dof + 1); 159424121865SAdrian Maldonado glob2loc[i] = dof; 159524121865SAdrian Maldonado } 159624121865SAdrian Maldonado 159724121865SAdrian Maldonado ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 159824121865SAdrian Maldonado #if 0 159924121865SAdrian Maldonado ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 160024121865SAdrian Maldonado #endif 160124121865SAdrian Maldonado PetscFunctionReturn(0); 160224121865SAdrian Maldonado } 160324121865SAdrian Maldonado 160401ad2aeeSHong Zhang #include <petsc/private/matimpl.h> 16051ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 16061ad426b7SHong Zhang { 16071ad426b7SHong Zhang PetscErrorCode ierr; 160824121865SAdrian Maldonado PetscMPIInt rank, size; 16091ad426b7SHong Zhang DM_Network *network = (DM_Network*) dm->data; 1610a4e85ca8SHong Zhang PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 1611840c2264SHong Zhang PetscInt cstart,ncols,j,e,v; 161224121865SAdrian Maldonado PetscBool ghost,ghost_vc,ghost2,isNest; 1613a4e85ca8SHong Zhang Mat Juser; 1614bfbc38dcSHong Zhang PetscSection sectionGlobal; 1615447d78afSSatish Balay PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 1616a4e85ca8SHong Zhang const PetscInt *edges,*cone; 16175cf7da58SHong Zhang MPI_Comm comm; 161824121865SAdrian Maldonado MatType mtype; 16195cf7da58SHong Zhang Vec vd_nz,vo_nz; 16205cf7da58SHong Zhang PetscInt *dnnz,*onnz; 16215cf7da58SHong Zhang PetscScalar *vdnz,*vonz; 16221ad426b7SHong Zhang 16231ad426b7SHong Zhang PetscFunctionBegin; 162424121865SAdrian Maldonado mtype = dm->mattype; 162524121865SAdrian Maldonado ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr); 162624121865SAdrian Maldonado 162724121865SAdrian Maldonado if (isNest) { 16280731d606SHong Zhang /* ierr = DMCreateMatrix_Network_Nest(); */ 162924121865SAdrian Maldonado PetscInt eDof, vDof; 163024121865SAdrian Maldonado Mat j11, j12, j21, j22, bA[2][2]; 163124121865SAdrian Maldonado ISLocalToGlobalMapping eISMap, vISMap; 163224121865SAdrian Maldonado 163324121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 163424121865SAdrian Maldonado ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 163524121865SAdrian Maldonado ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 163624121865SAdrian Maldonado 163724121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr); 163824121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr); 163924121865SAdrian Maldonado 164001ad2aeeSHong Zhang ierr = MatCreate(comm, &j11);CHKERRQ(ierr); 164124121865SAdrian Maldonado ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 164224121865SAdrian Maldonado ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr); 164324121865SAdrian Maldonado 164401ad2aeeSHong Zhang ierr = MatCreate(comm, &j12);CHKERRQ(ierr); 164524121865SAdrian Maldonado ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr); 164624121865SAdrian Maldonado ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr); 164724121865SAdrian Maldonado 164801ad2aeeSHong Zhang ierr = MatCreate(comm, &j21);CHKERRQ(ierr); 164924121865SAdrian Maldonado ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 165024121865SAdrian Maldonado ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr); 165124121865SAdrian Maldonado 165201ad2aeeSHong Zhang ierr = MatCreate(comm, &j22);CHKERRQ(ierr); 165324121865SAdrian Maldonado ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 165424121865SAdrian Maldonado ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr); 165524121865SAdrian Maldonado 16563f6a6bdaSHong Zhang bA[0][0] = j11; 16573f6a6bdaSHong Zhang bA[0][1] = j12; 16583f6a6bdaSHong Zhang bA[1][0] = j21; 16593f6a6bdaSHong Zhang bA[1][1] = j22; 166024121865SAdrian Maldonado 166124121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr); 166224121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr); 166324121865SAdrian Maldonado 166424121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr); 166524121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr); 166624121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr); 166724121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr); 166824121865SAdrian Maldonado 166924121865SAdrian Maldonado ierr = MatSetUp(j11);CHKERRQ(ierr); 167024121865SAdrian Maldonado ierr = MatSetUp(j12);CHKERRQ(ierr); 167124121865SAdrian Maldonado ierr = MatSetUp(j21);CHKERRQ(ierr); 167224121865SAdrian Maldonado ierr = MatSetUp(j22);CHKERRQ(ierr); 167324121865SAdrian Maldonado 167401ad2aeeSHong Zhang ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr); 167524121865SAdrian Maldonado ierr = MatSetUp(*J);CHKERRQ(ierr); 167624121865SAdrian Maldonado ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr); 167724121865SAdrian Maldonado ierr = MatDestroy(&j11);CHKERRQ(ierr); 167824121865SAdrian Maldonado ierr = MatDestroy(&j12);CHKERRQ(ierr); 167924121865SAdrian Maldonado ierr = MatDestroy(&j21);CHKERRQ(ierr); 168024121865SAdrian Maldonado ierr = MatDestroy(&j22);CHKERRQ(ierr); 168124121865SAdrian Maldonado 168224121865SAdrian Maldonado ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168324121865SAdrian Maldonado ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168424121865SAdrian Maldonado 168524121865SAdrian Maldonado /* Free structures */ 168624121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr); 168724121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr); 168824121865SAdrian Maldonado 168924121865SAdrian Maldonado PetscFunctionReturn(0); 169024121865SAdrian Maldonado } else if (!network->userEdgeJacobian && !network->userVertexJacobian) { 1691a4e85ca8SHong Zhang /* user does not provide Jacobian blocks */ 1692bfbc38dcSHong Zhang ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 1693bfbc38dcSHong Zhang ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 16941ad426b7SHong Zhang PetscFunctionReturn(0); 16951ad426b7SHong Zhang } 16961ad426b7SHong Zhang 1697bfbc38dcSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 16982a945128SHong Zhang ierr = DMGetDefaultGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 1699bfbc38dcSHong Zhang ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 1700bfbc38dcSHong Zhang ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 17012a945128SHong Zhang 17022a945128SHong Zhang ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 17032a945128SHong Zhang ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 170489898e50SHong Zhang 170589898e50SHong Zhang /* (1) Set matrix preallocation */ 170689898e50SHong Zhang /*------------------------------*/ 1707840c2264SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1708840c2264SHong Zhang ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr); 1709840c2264SHong Zhang ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr); 1710840c2264SHong Zhang ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr); 1711840c2264SHong Zhang ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr); 1712840c2264SHong Zhang ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr); 1713840c2264SHong Zhang 171489898e50SHong Zhang /* Set preallocation for edges */ 171589898e50SHong Zhang /*-----------------------------*/ 1716840c2264SHong Zhang ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 1717840c2264SHong Zhang 1718bdcb62a2SHong Zhang ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr); 1719840c2264SHong Zhang for (e=eStart; e<eEnd; e++) { 1720840c2264SHong Zhang /* Get row indices */ 1721840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 1722840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 1723840c2264SHong Zhang if (nrows) { 1724840c2264SHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 1725840c2264SHong Zhang 17265cf7da58SHong Zhang /* Set preallocation for conntected vertices */ 1727d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1728840c2264SHong Zhang for (v=0; v<2; v++) { 1729840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 1730840c2264SHong Zhang 17318675203cSHong Zhang if (network->Je) { 1732840c2264SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 17338675203cSHong Zhang } else Juser = NULL; 1734840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 17355cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1736840c2264SHong Zhang } 1737840c2264SHong Zhang 173889898e50SHong Zhang /* Set preallocation for edge self */ 1739840c2264SHong Zhang cstart = rstart; 17408675203cSHong Zhang if (network->Je) { 1741840c2264SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 17428675203cSHong Zhang } else Juser = NULL; 17435cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1744840c2264SHong Zhang } 1745840c2264SHong Zhang } 1746840c2264SHong Zhang 174789898e50SHong Zhang /* Set preallocation for vertices */ 174889898e50SHong Zhang /*--------------------------------*/ 1749840c2264SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 17508675203cSHong Zhang if (vEnd - vStart) vptr = network->Jvptr; 1751840c2264SHong Zhang 1752840c2264SHong Zhang for (v=vStart; v<vEnd; v++) { 1753840c2264SHong Zhang /* Get row indices */ 1754840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1755840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 1756840c2264SHong Zhang if (!nrows) continue; 1757840c2264SHong Zhang 1758bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1759bdcb62a2SHong Zhang if (ghost) { 1760bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1761bdcb62a2SHong Zhang } else { 1762bdcb62a2SHong Zhang rows_v = rows; 1763bdcb62a2SHong Zhang } 1764bdcb62a2SHong Zhang 1765bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1766840c2264SHong Zhang 1767840c2264SHong Zhang /* Get supporting edges and connected vertices */ 1768840c2264SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1769840c2264SHong Zhang 1770840c2264SHong Zhang for (e=0; e<nedges; e++) { 1771840c2264SHong Zhang /* Supporting edges */ 1772840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1773840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1774840c2264SHong Zhang 17758675203cSHong Zhang if (network->Jv) { 1776840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 17778675203cSHong Zhang } else Juser = NULL; 1778bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1779840c2264SHong Zhang 1780840c2264SHong Zhang /* Connected vertices */ 1781d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 1782840c2264SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 1783840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr); 1784840c2264SHong Zhang 1785840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1786840c2264SHong Zhang 17878675203cSHong Zhang if (network->Jv) { 1788840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 17898675203cSHong Zhang } else Juser = NULL; 1790e102a522SHong Zhang if (ghost_vc||ghost) { 1791e102a522SHong Zhang ghost2 = PETSC_TRUE; 1792e102a522SHong Zhang } else { 1793e102a522SHong Zhang ghost2 = PETSC_FALSE; 1794e102a522SHong Zhang } 1795e102a522SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr); 1796840c2264SHong Zhang } 1797840c2264SHong Zhang 179889898e50SHong Zhang /* Set preallocation for vertex self */ 1799840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1800840c2264SHong Zhang if (!ghost) { 1801840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 18028675203cSHong Zhang if (network->Jv) { 1803840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 18048675203cSHong Zhang } else Juser = NULL; 1805bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1806840c2264SHong Zhang } 1807bdcb62a2SHong Zhang if (ghost) { 1808bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1809bdcb62a2SHong Zhang } 1810840c2264SHong Zhang } 1811840c2264SHong Zhang 1812840c2264SHong Zhang ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr); 1813840c2264SHong Zhang ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr); 18145cf7da58SHong Zhang 18155cf7da58SHong Zhang ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr); 18165cf7da58SHong Zhang 18175cf7da58SHong Zhang ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr); 1818840c2264SHong Zhang ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr); 1819840c2264SHong Zhang 1820840c2264SHong Zhang ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr); 1821840c2264SHong Zhang ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr); 1822840c2264SHong Zhang for (j=0; j<localSize; j++) { 1823e102a522SHong Zhang dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 1824e102a522SHong Zhang onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 1825840c2264SHong Zhang } 1826840c2264SHong Zhang ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr); 1827840c2264SHong Zhang ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr); 1828840c2264SHong Zhang ierr = VecDestroy(&vd_nz);CHKERRQ(ierr); 1829840c2264SHong Zhang ierr = VecDestroy(&vo_nz);CHKERRQ(ierr); 1830840c2264SHong Zhang 18315cf7da58SHong Zhang ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr); 18325cf7da58SHong Zhang ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr); 18335cf7da58SHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 18345cf7da58SHong Zhang 18355cf7da58SHong Zhang ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr); 18365cf7da58SHong Zhang 183789898e50SHong Zhang /* (2) Set matrix entries for edges */ 183889898e50SHong Zhang /*----------------------------------*/ 18391ad426b7SHong Zhang for (e=eStart; e<eEnd; e++) { 1840bfbc38dcSHong Zhang /* Get row indices */ 18411ad426b7SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 184217df6e9eSHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 18434b976069SHong Zhang if (nrows) { 184417df6e9eSHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 18451ad426b7SHong Zhang 1846bfbc38dcSHong Zhang /* Set matrix entries for conntected vertices */ 1847d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1848bfbc38dcSHong Zhang for (v=0; v<2; v++) { 1849bfbc38dcSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 1850883e35e8SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 18513e97b6e8SHong Zhang 18528675203cSHong Zhang if (network->Je) { 1853a4e85ca8SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 18548675203cSHong Zhang } else Juser = NULL; 1855a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1856bfbc38dcSHong Zhang } 185717df6e9eSHong Zhang 1858bfbc38dcSHong Zhang /* Set matrix entries for edge self */ 18593e97b6e8SHong Zhang cstart = rstart; 18608675203cSHong Zhang if (network->Je) { 1861a4e85ca8SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 18628675203cSHong Zhang } else Juser = NULL; 1863a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 18641ad426b7SHong Zhang } 18654b976069SHong Zhang } 18661ad426b7SHong Zhang 1867bfbc38dcSHong Zhang /* Set matrix entries for vertices */ 186883b2e829SHong Zhang /*---------------------------------*/ 18691ad426b7SHong Zhang for (v=vStart; v<vEnd; v++) { 1870bfbc38dcSHong Zhang /* Get row indices */ 1871596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1872596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 18734b976069SHong Zhang if (!nrows) continue; 1874596e729fSHong Zhang 1875bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1876bdcb62a2SHong Zhang if (ghost) { 1877bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1878bdcb62a2SHong Zhang } else { 1879bdcb62a2SHong Zhang rows_v = rows; 1880bdcb62a2SHong Zhang } 1881bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1882596e729fSHong Zhang 1883bfbc38dcSHong Zhang /* Get supporting edges and connected vertices */ 1884596e729fSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1885596e729fSHong Zhang 1886596e729fSHong Zhang for (e=0; e<nedges; e++) { 1887bfbc38dcSHong Zhang /* Supporting edges */ 1888596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1889596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1890596e729fSHong Zhang 18918675203cSHong Zhang if (network->Jv) { 1892a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 18938675203cSHong Zhang } else Juser = NULL; 1894bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1895596e729fSHong Zhang 1896bfbc38dcSHong Zhang /* Connected vertices */ 1897d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 18982a945128SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 18992a945128SHong Zhang 190044aca652SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 190144aca652SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1902a4e85ca8SHong Zhang 19038675203cSHong Zhang if (network->Jv) { 1904a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 19058675203cSHong Zhang } else Juser = NULL; 1906bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1907596e729fSHong Zhang } 1908596e729fSHong Zhang 1909bfbc38dcSHong Zhang /* Set matrix entries for vertex self */ 19101ad426b7SHong Zhang if (!ghost) { 1911596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 19128675203cSHong Zhang if (network->Jv) { 1913a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 19148675203cSHong Zhang } else Juser = NULL; 1915bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr); 1916bdcb62a2SHong Zhang } 1917bdcb62a2SHong Zhang if (ghost) { 1918bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1919bdcb62a2SHong Zhang } 19201ad426b7SHong Zhang } 1921a4e85ca8SHong Zhang ierr = PetscFree(rows);CHKERRQ(ierr); 1922bdcb62a2SHong Zhang 19231ad426b7SHong Zhang ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19241ad426b7SHong Zhang ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1925dd6f46cdSHong Zhang 19265f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 19275f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19285f2c45f1SShri Abhyankar } 19295f2c45f1SShri Abhyankar 19305f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 19315f2c45f1SShri Abhyankar { 19325f2c45f1SShri Abhyankar PetscErrorCode ierr; 19335f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19342727e31bSShri Abhyankar PetscInt j; 19355f2c45f1SShri Abhyankar 19365f2c45f1SShri Abhyankar PetscFunctionBegin; 19378415c774SShri Abhyankar if (--network->refct > 0) PetscFunctionReturn(0); 193883b2e829SHong Zhang if (network->Je) { 193983b2e829SHong Zhang ierr = PetscFree(network->Je);CHKERRQ(ierr); 194083b2e829SHong Zhang } 194183b2e829SHong Zhang if (network->Jv) { 1942883e35e8SHong Zhang ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 194383b2e829SHong Zhang ierr = PetscFree(network->Jv);CHKERRQ(ierr); 19441ad426b7SHong Zhang } 194513c2a604SAdrian Maldonado 194613c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr); 194713c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr); 194813c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr); 194913c2a604SAdrian Maldonado if (network->vertex.sf) { 195013c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr); 195113c2a604SAdrian Maldonado } 195213c2a604SAdrian Maldonado /* edge */ 195313c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr); 195413c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr); 195513c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr); 195613c2a604SAdrian Maldonado if (network->edge.sf) { 195713c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr); 195813c2a604SAdrian Maldonado } 19595f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 19605f2c45f1SShri Abhyankar network->edges = NULL; 19615f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 19625f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 196383b2e829SHong Zhang 19642727e31bSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 19652727e31bSShri Abhyankar ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr); 19662727e31bSShri Abhyankar ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr); 19672727e31bSShri Abhyankar } 1968e2aaf10cSShri Abhyankar ierr = PetscFree(network->subnet);CHKERRQ(ierr); 19695f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 19705f2c45f1SShri Abhyankar ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 19715f2c45f1SShri Abhyankar ierr = PetscFree(network->header);CHKERRQ(ierr); 19725f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 19735f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19745f2c45f1SShri Abhyankar } 19755f2c45f1SShri Abhyankar 19765f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 19775f2c45f1SShri Abhyankar { 19785f2c45f1SShri Abhyankar PetscErrorCode ierr; 19795f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19805f2c45f1SShri Abhyankar 19815f2c45f1SShri Abhyankar PetscFunctionBegin; 19825f2c45f1SShri Abhyankar ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 19835f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19845f2c45f1SShri Abhyankar } 19855f2c45f1SShri Abhyankar 19865f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 19875f2c45f1SShri Abhyankar { 19885f2c45f1SShri Abhyankar PetscErrorCode ierr; 19895f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19905f2c45f1SShri Abhyankar 19915f2c45f1SShri Abhyankar PetscFunctionBegin; 19925f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 19935f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19945f2c45f1SShri Abhyankar } 19955f2c45f1SShri Abhyankar 19965f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 19975f2c45f1SShri Abhyankar { 19985f2c45f1SShri Abhyankar PetscErrorCode ierr; 19995f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 20005f2c45f1SShri Abhyankar 20015f2c45f1SShri Abhyankar PetscFunctionBegin; 20025f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 20035f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20045f2c45f1SShri Abhyankar } 20055f2c45f1SShri Abhyankar 20065f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 20075f2c45f1SShri Abhyankar { 20085f2c45f1SShri Abhyankar PetscErrorCode ierr; 20095f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 20105f2c45f1SShri Abhyankar 20115f2c45f1SShri Abhyankar PetscFunctionBegin; 20125f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 20135f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20145f2c45f1SShri Abhyankar } 20155f2c45f1SShri Abhyankar 20165f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 20175f2c45f1SShri Abhyankar { 20185f2c45f1SShri Abhyankar PetscErrorCode ierr; 20195f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 20205f2c45f1SShri Abhyankar 20215f2c45f1SShri Abhyankar PetscFunctionBegin; 20225f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 20235f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20245f2c45f1SShri Abhyankar } 2025