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 35*f025b11dSHong Zhang . NsubnetCouple - number of coupling subnetworks 36*f025b11dSHong Zhang . nV - number of local vertices for each subnetwork and coupling subnetwork 37*f025b11dSHong Zhang . nE - number of local edges for each subnetwork and coupling subnetwork 38*f025b11dSHong Zhang . NV - number of global vertices (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork 39*f025b11dSHong Zhang - NE - number of global edges (or PETSC_DETERMINE) for each subnetwork and coupling subnetwork 405f2c45f1SShri Abhyankar 415f2c45f1SShri Abhyankar Notes 425f2c45f1SShri Abhyankar If one processor calls this with NV (NE) of PETSC_DECIDE then all processors must, otherwise the prgram will hang. 435f2c45f1SShri Abhyankar 44*f025b11dSHong Zhang You cannot change the sizes once they have been set. nV, nE, NV and NE are arrays of length Nsubnet+NsubnetCouple. 455f2c45f1SShri Abhyankar 461b266c99SBarry Smith Level: intermediate 471b266c99SBarry Smith 481b266c99SBarry Smith .seealso: DMNetworkCreate() 495f2c45f1SShri Abhyankar @*/ 50*f025b11dSHong Zhang PetscErrorCode DMNetworkSetSizes(DM dm,PetscInt Nsubnet,PetscInt NsubnetCouple,PetscInt nV[], PetscInt nE[], PetscInt NV[], PetscInt NE[]) 515f2c45f1SShri Abhyankar { 525f2c45f1SShri Abhyankar PetscErrorCode ierr; 535f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 54e2aaf10cSShri Abhyankar PetscInt a[2],b[2],i; 555f2c45f1SShri Abhyankar 565f2c45f1SShri Abhyankar PetscFunctionBegin; 575f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 58e2aaf10cSShri Abhyankar if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet); 597765340cSHong Zhang if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple); 607765340cSHong Zhang 617765340cSHong Zhang if (Nsubnet > 0) PetscValidLogicalCollectiveInt(dm,Nsubnet,2); 627765340cSHong Zhang if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network"); 637765340cSHong Zhang 647765340cSHong Zhang network->nsubnet = Nsubnet + NsubnetCouple; 657765340cSHong Zhang Nsubnet += NsubnetCouple; 667765340cSHong Zhang network->ncsubnet = NsubnetCouple; 677765340cSHong Zhang ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr); 687765340cSHong Zhang for(i=0; i < Nsubnet; i++) { 69*f025b11dSHong Zhang if (NV) { 707765340cSHong Zhang if (NV[i] > 0) PetscValidLogicalCollectiveInt(dm,NV[i],6); 717765340cSHong 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]); 72*f025b11dSHong Zhang } 73*f025b11dSHong Zhang if (NE) { 74*f025b11dSHong Zhang if (NE[i] > 0) PetscValidLogicalCollectiveInt(dm,NE[i],7); 757765340cSHong 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]); 76*f025b11dSHong Zhang } 77*f025b11dSHong Zhang 787765340cSHong Zhang a[0] = nV[i]; a[1] = nE[i]; 797765340cSHong Zhang ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 807765340cSHong Zhang network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1]; 817765340cSHong Zhang 827765340cSHong Zhang network->subnet[i].id = i; 837765340cSHong Zhang 847765340cSHong Zhang network->subnet[i].nvtx = nV[i]; 857765340cSHong Zhang network->subnet[i].vStart = network->nVertices; 867765340cSHong Zhang network->subnet[i].vEnd = network->subnet[i].vStart + network->subnet[i].Nvtx; 877765340cSHong Zhang network->nVertices += network->subnet[i].nvtx; 887765340cSHong Zhang network->NVertices += network->subnet[i].Nvtx; 897765340cSHong Zhang 907765340cSHong Zhang network->subnet[i].nedge = nE[i]; 917765340cSHong Zhang network->subnet[i].eStart = network->nEdges; 927765340cSHong Zhang network->subnet[i].eEnd = network->subnet[i].eStart + network->subnet[i].Nedge; 937765340cSHong Zhang network->nEdges += network->subnet[i].nedge; 947765340cSHong Zhang network->NEdges += network->subnet[i].Nedge; 957765340cSHong Zhang } 967765340cSHong Zhang PetscFunctionReturn(0); 977765340cSHong Zhang } 987765340cSHong Zhang 995f2c45f1SShri Abhyankar /*@ 1005f2c45f1SShri Abhyankar DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 1015f2c45f1SShri Abhyankar 1025f2c45f1SShri Abhyankar Logically collective on DM 1035f2c45f1SShri Abhyankar 1045f2c45f1SShri Abhyankar Input Parameters: 105e2aaf10cSShri Abhyankar . edges - list of edges for each subnetwork 1065f2c45f1SShri Abhyankar 1075f2c45f1SShri Abhyankar Notes: 1085f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 1095f2c45f1SShri Abhyankar not be destroyed before the call to DMNetworkLayoutSetUp 1105f2c45f1SShri Abhyankar 1115f2c45f1SShri Abhyankar Level: intermediate 1125f2c45f1SShri Abhyankar 1135f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes 1145f2c45f1SShri Abhyankar @*/ 115e2aaf10cSShri Abhyankar PetscErrorCode DMNetworkSetEdgeList(DM dm, int *edgelist[]) 1165f2c45f1SShri Abhyankar { 1175f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 118e2aaf10cSShri Abhyankar PetscInt i; 1195f2c45f1SShri Abhyankar 1205f2c45f1SShri Abhyankar PetscFunctionBegin; 121e2aaf10cSShri Abhyankar for(i=0; i < network->nsubnet; i++) { 122e2aaf10cSShri Abhyankar network->subnet[i].edgelist = edgelist[i]; 123e2aaf10cSShri Abhyankar } 1245f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1255f2c45f1SShri Abhyankar } 126*f025b11dSHong Zhang #if 0 127*f025b11dSHong Zhang PetscErrorCode DMNetworkSetEdgeListCoupled(DM dm,int *edgelist[],int *edgelistCouple[]) 128*f025b11dSHong Zhang { 129*f025b11dSHong Zhang DM_Network *network = (DM_Network*) dm->data; 130*f025b11dSHong Zhang PetscInt i; 131*f025b11dSHong Zhang 132*f025b11dSHong Zhang PetscFunctionBegin; 133*f025b11dSHong Zhang for(i=0; i < network->nsubnet; i++) { 134*f025b11dSHong Zhang network->subnet[i].edgelist = edgelist[i]; 135*f025b11dSHong Zhang } 136*f025b11dSHong Zhang PetscFunctionReturn(0); 137*f025b11dSHong Zhang } 138*f025b11dSHong Zhang #endif 1395f2c45f1SShri Abhyankar 1405f2c45f1SShri Abhyankar /*@ 1415f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 1425f2c45f1SShri Abhyankar 1435f2c45f1SShri Abhyankar Collective on DM 1445f2c45f1SShri Abhyankar 1455f2c45f1SShri Abhyankar Input Parameters 1465f2c45f1SShri Abhyankar . DM - the dmnetwork object 1475f2c45f1SShri Abhyankar 1485f2c45f1SShri Abhyankar Notes: 1495f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 1505f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 1515f2c45f1SShri Abhyankar 1525f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 1535f2c45f1SShri Abhyankar 1545f2c45f1SShri Abhyankar Level: intermediate 1555f2c45f1SShri Abhyankar 1565f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 1575f2c45f1SShri Abhyankar @*/ 1585f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 1595f2c45f1SShri Abhyankar { 1605f2c45f1SShri Abhyankar PetscErrorCode ierr; 1615f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 1625f2c45f1SShri Abhyankar PetscInt dim = 1; /* One dimensional network */ 1635f2c45f1SShri Abhyankar PetscInt numCorners=2; 1645f2c45f1SShri Abhyankar PetscInt spacedim=2; 1655f2c45f1SShri Abhyankar double *vertexcoords=NULL; 166e2aaf10cSShri Abhyankar PetscInt i,j; 1675f2c45f1SShri Abhyankar PetscInt ndata; 168e2aaf10cSShri Abhyankar PetscInt ctr=0; 1695f2c45f1SShri Abhyankar 1705f2c45f1SShri Abhyankar PetscFunctionBegin; 1716fefedf4SHong Zhang if (network->nVertices) { 1726fefedf4SHong Zhang ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr); 1735f2c45f1SShri Abhyankar } 174e2aaf10cSShri Abhyankar 175e2aaf10cSShri Abhyankar /* Create the edgelist for the network by concatenating edgelists of the subnetworks */ 176e2aaf10cSShri Abhyankar ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr); 177e2aaf10cSShri Abhyankar for(i=0; i < network->nsubnet; i++) { 178e2aaf10cSShri Abhyankar for(j = 0; j < network->subnet[i].nedge; j++) { 179e2aaf10cSShri Abhyankar network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j]; 180e2aaf10cSShri Abhyankar network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1]; 181e2aaf10cSShri Abhyankar ctr++; 182e2aaf10cSShri Abhyankar } 183e2aaf10cSShri Abhyankar } 184e2aaf10cSShri Abhyankar 18561de3474SHong Zhang #if 0 186e2aaf10cSShri Abhyankar for(i=0; i < network->nEdges; i++) { 187e2aaf10cSShri Abhyankar ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr); 188e2aaf10cSShri Abhyankar } 189e2aaf10cSShri Abhyankar #endif 190e2aaf10cSShri Abhyankar 1916fefedf4SHong Zhang ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 1926fefedf4SHong Zhang if (network->nVertices) { 1935f2c45f1SShri Abhyankar ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 1945f2c45f1SShri Abhyankar } 195e2aaf10cSShri Abhyankar ierr = PetscFree(network->edges);CHKERRQ(ierr); 196e2aaf10cSShri Abhyankar 1975f2c45f1SShri Abhyankar ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 1985f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 1995f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 2005f2c45f1SShri Abhyankar 2015f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 2025f2c45f1SShri Abhyankar ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 2035f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 2045f2c45f1SShri Abhyankar ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 2055f2c45f1SShri Abhyankar 2062727e31bSShri Abhyankar /* Create vertices and edges array for the subnetworks */ 2072727e31bSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 2082727e31bSShri Abhyankar ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr); 2092727e31bSShri Abhyankar ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr); 2102727e31bSShri Abhyankar /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 2112727e31bSShri Abhyankar These get updated when the vertices and edges are added. */ 2122727e31bSShri Abhyankar network->subnet[j].nvtx = network->subnet[j].nedge = 0; 2132727e31bSShri Abhyankar } 2142727e31bSShri Abhyankar 2155f2c45f1SShri Abhyankar network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 2166caa05f4SBarry Smith ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 217e2aaf10cSShri Abhyankar for(i=network->eStart; i < network->eEnd; i++) { 218e2aaf10cSShri Abhyankar network->header[i].index = i; /* Global edge number */ 219e2aaf10cSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 220e2aaf10cSShri Abhyankar if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) { 221e2aaf10cSShri Abhyankar network->header[i].subnetid = j; /* Subnetwork id */ 2222727e31bSShri Abhyankar network->subnet[j].edges[network->subnet[j].nedge++] = i; 223e2aaf10cSShri Abhyankar break; 224e2aaf10cSShri Abhyankar } 2257b6afd5bSHong Zhang } 2265f2c45f1SShri Abhyankar network->header[i].ndata = 0; 2275f2c45f1SShri Abhyankar ndata = network->header[i].ndata; 2285f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 2295f2c45f1SShri Abhyankar network->header[i].offset[ndata] = 0; 2305f2c45f1SShri Abhyankar } 231e2aaf10cSShri Abhyankar 232e2aaf10cSShri Abhyankar for(i=network->vStart; i < network->vEnd; i++) { 233e2aaf10cSShri Abhyankar network->header[i].index = i - network->vStart; 234e2aaf10cSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 235e2aaf10cSShri Abhyankar if((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) { 236e2aaf10cSShri Abhyankar network->header[i].subnetid = j; 2372727e31bSShri Abhyankar network->subnet[j].vertices[network->subnet[j].nvtx++] = i; 238e2aaf10cSShri Abhyankar break; 239e2aaf10cSShri Abhyankar } 240e2aaf10cSShri Abhyankar } 241e2aaf10cSShri Abhyankar network->header[i].ndata = 0; 242e2aaf10cSShri Abhyankar ndata = network->header[i].ndata; 243e2aaf10cSShri Abhyankar ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 244e2aaf10cSShri Abhyankar network->header[i].offset[ndata] = 0; 245e2aaf10cSShri Abhyankar } 246e2aaf10cSShri Abhyankar 247854ce69bSBarry Smith ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 2486500d4abSHong Zhang PetscFunctionReturn(0); 2496500d4abSHong Zhang } 250e2aaf10cSShri Abhyankar 2516500d4abSHong Zhang PetscErrorCode DMNetworkLayoutSetUpCoupled(DM dm) 2526500d4abSHong Zhang { 2536500d4abSHong Zhang PetscErrorCode ierr; 2546500d4abSHong Zhang DM_Network *network = (DM_Network*)dm->data; 2556500d4abSHong Zhang PetscInt dim = 1; /* One dimensional network */ 256991cf414SHong Zhang PetscInt numCorners=2,spacedim=2; 2576500d4abSHong Zhang double *vertexcoords=NULL; 2587765340cSHong Zhang PetscInt i,j,ndata,ctr=0,nsubnet; 259991cf414SHong Zhang PetscInt *edgelist_couple=NULL,k,netid,vid; 2606500d4abSHong Zhang 2616500d4abSHong Zhang PetscFunctionBegin; 2626500d4abSHong Zhang if (network->nVertices) { 2636500d4abSHong Zhang ierr = PetscCalloc1(numCorners*network->nVertices,&vertexcoords);CHKERRQ(ierr); 2646500d4abSHong Zhang } 2656500d4abSHong Zhang 2666500d4abSHong Zhang /* Create the edgelist for the network by concatenating edgelists of the subnetworks */ 2677765340cSHong Zhang nsubnet = network->nsubnet - network->ncsubnet; 2686500d4abSHong Zhang ierr = PetscCalloc1(2*network->nEdges,&network->edges);CHKERRQ(ierr); 2697765340cSHong Zhang for (i=0; i < nsubnet; i++) { 2706500d4abSHong Zhang for (j = 0; j < network->subnet[i].nedge; j++) { 2716500d4abSHong Zhang network->edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j]; 2726500d4abSHong Zhang network->edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1]; 2736500d4abSHong Zhang ctr++; 2746500d4abSHong Zhang } 2756500d4abSHong Zhang } 2767765340cSHong Zhang 2777765340cSHong Zhang i = nsubnet; /* netid of coupling subnet */ 2787765340cSHong Zhang nsubnet = network->nsubnet; 2797765340cSHong Zhang while (i < nsubnet) { 280991cf414SHong Zhang edgelist_couple = network->subnet[i].edgelist; 2816500d4abSHong Zhang k = 0; 2826500d4abSHong Zhang for (j = 0; j < network->subnet[i].nedge; j++) { 2836500d4abSHong Zhang netid = edgelist_couple[k]; vid = edgelist_couple[k+1]; 2846500d4abSHong Zhang network->edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2; 2856500d4abSHong Zhang 2866500d4abSHong Zhang netid = edgelist_couple[k]; vid = edgelist_couple[k+1]; 287991cf414SHong Zhang network->edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2; 2886500d4abSHong Zhang ctr++; 2896500d4abSHong Zhang } 2907765340cSHong Zhang i++; 2917765340cSHong Zhang } 2926500d4abSHong Zhang 29334bd1bf5SHong Zhang #if 0 2946500d4abSHong Zhang for(i=0; i < network->nEdges; i++) { 2956500d4abSHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",network->edges[2*i],network->edges[2*i+1]);CHKERRQ(ierr); 2966500d4abSHong Zhang printf("\n"); 2976500d4abSHong Zhang } 2986500d4abSHong Zhang #endif 2996500d4abSHong Zhang 3006500d4abSHong Zhang ierr = DMPlexCreateFromCellList(PetscObjectComm((PetscObject)dm),dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,network->edges,spacedim,vertexcoords,&network->plex);CHKERRQ(ierr); 3016500d4abSHong Zhang if (network->nVertices) { 3026500d4abSHong Zhang ierr = PetscFree(vertexcoords);CHKERRQ(ierr); 3036500d4abSHong Zhang } 3046500d4abSHong Zhang ierr = PetscFree(network->edges);CHKERRQ(ierr); 3056500d4abSHong Zhang 3066500d4abSHong Zhang ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 3076500d4abSHong Zhang ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 3086500d4abSHong Zhang ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 3096500d4abSHong Zhang 3106500d4abSHong Zhang ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DataSection);CHKERRQ(ierr); 3116500d4abSHong Zhang ierr = PetscSectionCreate(PetscObjectComm((PetscObject)dm),&network->DofSection);CHKERRQ(ierr); 3126500d4abSHong Zhang ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 3136500d4abSHong Zhang ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 3146500d4abSHong Zhang 3156500d4abSHong Zhang /* Create vertices and edges array for the subnetworks */ 3166500d4abSHong Zhang for (j=0; j < network->nsubnet; j++) { 3176500d4abSHong Zhang ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr); 3186500d4abSHong Zhang ierr = PetscCalloc1(network->subnet[j].nvtx,&network->subnet[j].vertices);CHKERRQ(ierr); 3196500d4abSHong Zhang /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 3206500d4abSHong Zhang These get updated when the vertices and edges are added. */ 3216500d4abSHong Zhang network->subnet[j].nvtx = network->subnet[j].nedge = 0; 3226500d4abSHong Zhang } 3236500d4abSHong Zhang 3246500d4abSHong Zhang network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 3256500d4abSHong Zhang ierr = PetscCalloc1(network->pEnd-network->pStart,&network->header);CHKERRQ(ierr); 3266500d4abSHong Zhang for (i=network->eStart; i < network->eEnd; i++) { 3276500d4abSHong Zhang network->header[i].index = i; /* Global edge number */ 3286500d4abSHong Zhang for (j=0; j < network->nsubnet; j++) { 3296500d4abSHong Zhang if((network->subnet[j].eStart <= i) && (i < network->subnet[j].eEnd)) { 3306500d4abSHong Zhang network->header[i].subnetid = j; /* Subnetwork id */ 3316500d4abSHong Zhang network->subnet[j].edges[network->subnet[j].nedge++] = i; 3326500d4abSHong Zhang break; 3336500d4abSHong Zhang } 3346500d4abSHong Zhang } 3356500d4abSHong Zhang network->header[i].ndata = 0; 3366500d4abSHong Zhang ndata = network->header[i].ndata; 3376500d4abSHong Zhang ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 3386500d4abSHong Zhang network->header[i].offset[ndata] = 0; 3396500d4abSHong Zhang } 3406500d4abSHong Zhang 3416500d4abSHong Zhang for(i=network->vStart; i < network->vEnd; i++) { 3426500d4abSHong Zhang network->header[i].index = i - network->vStart; 3436500d4abSHong Zhang for (j=0; j < network->nsubnet; j++) { 3446500d4abSHong Zhang if ((network->subnet[j].vStart <= i-network->vStart) && (i-network->vStart < network->subnet[j].vEnd)) { 3456500d4abSHong Zhang network->header[i].subnetid = j; 3466500d4abSHong Zhang network->subnet[j].vertices[network->subnet[j].nvtx++] = i; 3476500d4abSHong Zhang break; 3486500d4abSHong Zhang } 3496500d4abSHong Zhang } 3506500d4abSHong Zhang network->header[i].ndata = 0; 3516500d4abSHong Zhang ndata = network->header[i].ndata; 3526500d4abSHong Zhang ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 3536500d4abSHong Zhang network->header[i].offset[ndata] = 0; 3546500d4abSHong Zhang } 3556500d4abSHong Zhang 3566500d4abSHong Zhang ierr = PetscMalloc1(network->pEnd-network->pStart,&network->cvalue);CHKERRQ(ierr); 3575f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3585f2c45f1SShri Abhyankar } 3595f2c45f1SShri Abhyankar 36094ef8ddeSSatish Balay /*@C 3612727e31bSShri Abhyankar DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork 3622727e31bSShri Abhyankar 3632727e31bSShri Abhyankar Input Parameters 3642727e31bSShri Abhyankar + dm - the number object 3652727e31bSShri Abhyankar - id - the ID (integer) of the subnetwork 3662727e31bSShri Abhyankar 3672727e31bSShri Abhyankar Output Parameters 3682727e31bSShri Abhyankar + nv - number of vertices (local) 3692727e31bSShri Abhyankar . ne - number of edges (local) 3702727e31bSShri Abhyankar . vtx - local vertices for this subnetwork 3712727e31bSShri Abhyankar . edge - local edges for this subnetwork 3722727e31bSShri Abhyankar 3732727e31bSShri Abhyankar Notes: 3742727e31bSShri Abhyankar Cannot call this routine before DMNetworkLayoutSetup() 3752727e31bSShri Abhyankar 3762727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 3772727e31bSShri Abhyankar @*/ 3782727e31bSShri Abhyankar PetscErrorCode DMNetworkGetSubnetworkInfo(DM netdm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge) 3792727e31bSShri Abhyankar { 3802727e31bSShri Abhyankar DM_Network *network = (DM_Network*) netdm->data; 3812727e31bSShri Abhyankar 3822727e31bSShri Abhyankar PetscFunctionBegin; 3832727e31bSShri Abhyankar *nv = network->subnet[id].nvtx; 3842727e31bSShri Abhyankar *ne = network->subnet[id].nedge; 3852727e31bSShri Abhyankar *vtx = network->subnet[id].vertices; 3862727e31bSShri Abhyankar *edge = network->subnet[id].edges; 3872727e31bSShri Abhyankar PetscFunctionReturn(0); 3882727e31bSShri Abhyankar } 3892727e31bSShri Abhyankar 3902727e31bSShri Abhyankar /*@C 3915f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 3925f2c45f1SShri Abhyankar 3935f2c45f1SShri Abhyankar Logically collective on DM 3945f2c45f1SShri Abhyankar 3955f2c45f1SShri Abhyankar Input Parameters 3965f2c45f1SShri Abhyankar + dm - the network object 3975f2c45f1SShri Abhyankar . name - the component name 3985f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 3995f2c45f1SShri Abhyankar 4005f2c45f1SShri Abhyankar Output Parameters 4015f2c45f1SShri Abhyankar . key - an integer key that defines the component 4025f2c45f1SShri Abhyankar 4035f2c45f1SShri Abhyankar Notes 4045f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 4055f2c45f1SShri Abhyankar 4065f2c45f1SShri Abhyankar Level: intermediate 4075f2c45f1SShri Abhyankar 4085f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 4095f2c45f1SShri Abhyankar @*/ 4105f2c45f1SShri Abhyankar PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,PetscInt size,PetscInt *key) 4115f2c45f1SShri Abhyankar { 4125f2c45f1SShri Abhyankar PetscErrorCode ierr; 4135f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 4145f2c45f1SShri Abhyankar DMNetworkComponent *component=&network->component[network->ncomponent]; 4155f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 4165f2c45f1SShri Abhyankar PetscInt i; 4175f2c45f1SShri Abhyankar 4185f2c45f1SShri Abhyankar PetscFunctionBegin; 4195f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 4205f2c45f1SShri Abhyankar ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 4215f2c45f1SShri Abhyankar if (flg) { 4225f2c45f1SShri Abhyankar *key = i; 4235f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4245f2c45f1SShri Abhyankar } 4256d64e262SShri Abhyankar } 4266d64e262SShri Abhyankar if(network->ncomponent == MAX_COMPONENTS) { 4276d64e262SShri Abhyankar SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS); 4285f2c45f1SShri Abhyankar } 4295f2c45f1SShri Abhyankar 4305f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 4315f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 4325f2c45f1SShri Abhyankar *key = network->ncomponent; 4335f2c45f1SShri Abhyankar network->ncomponent++; 4345f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4355f2c45f1SShri Abhyankar } 4365f2c45f1SShri Abhyankar 4375f2c45f1SShri Abhyankar /*@ 4385f2c45f1SShri Abhyankar DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 4395f2c45f1SShri Abhyankar 4405f2c45f1SShri Abhyankar Not Collective 4415f2c45f1SShri Abhyankar 4425f2c45f1SShri Abhyankar Input Parameters: 4435f2c45f1SShri Abhyankar + dm - The DMNetwork object 4445f2c45f1SShri Abhyankar 4455f2c45f1SShri Abhyankar Output Paramters: 4465f2c45f1SShri Abhyankar + vStart - The first vertex point 4475f2c45f1SShri Abhyankar - vEnd - One beyond the last vertex point 4485f2c45f1SShri Abhyankar 4495f2c45f1SShri Abhyankar Level: intermediate 4505f2c45f1SShri Abhyankar 4515f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange 4525f2c45f1SShri Abhyankar @*/ 4535f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 4545f2c45f1SShri Abhyankar { 4555f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4565f2c45f1SShri Abhyankar 4575f2c45f1SShri Abhyankar PetscFunctionBegin; 4585f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 4595f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 4605f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4615f2c45f1SShri Abhyankar } 4625f2c45f1SShri Abhyankar 4635f2c45f1SShri Abhyankar /*@ 4645f2c45f1SShri Abhyankar DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 4655f2c45f1SShri Abhyankar 4665f2c45f1SShri Abhyankar Not Collective 4675f2c45f1SShri Abhyankar 4685f2c45f1SShri Abhyankar Input Parameters: 4695f2c45f1SShri Abhyankar + dm - The DMNetwork object 4705f2c45f1SShri Abhyankar 4715f2c45f1SShri Abhyankar Output Paramters: 4725f2c45f1SShri Abhyankar + eStart - The first edge point 4735f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 4745f2c45f1SShri Abhyankar 4755f2c45f1SShri Abhyankar Level: intermediate 4765f2c45f1SShri Abhyankar 4775f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange 4785f2c45f1SShri Abhyankar @*/ 4795f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 4805f2c45f1SShri Abhyankar { 4815f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4825f2c45f1SShri Abhyankar 4835f2c45f1SShri Abhyankar PetscFunctionBegin; 4845f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 4855f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 4865f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4875f2c45f1SShri Abhyankar } 4885f2c45f1SShri Abhyankar 4897b6afd5bSHong Zhang /*@ 490e85e6aecSHong Zhang DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge. 4917b6afd5bSHong Zhang 4927b6afd5bSHong Zhang Not Collective 4937b6afd5bSHong Zhang 4947b6afd5bSHong Zhang Input Parameters: 4957b6afd5bSHong Zhang + dm - DMNetwork object 496e85e6aecSHong Zhang - p - edge point 4977b6afd5bSHong Zhang 4987b6afd5bSHong Zhang Output Paramters: 499e85e6aecSHong Zhang . index - user global numbering for the edge 5007b6afd5bSHong Zhang 5017b6afd5bSHong Zhang Level: intermediate 5027b6afd5bSHong Zhang 503e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex 5047b6afd5bSHong Zhang @*/ 505e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index) 5067b6afd5bSHong Zhang { 5077b6afd5bSHong Zhang PetscErrorCode ierr; 5087b6afd5bSHong Zhang DM_Network *network = (DM_Network*)dm->data; 5097b6afd5bSHong Zhang PetscInt offsetp; 5107b6afd5bSHong Zhang DMNetworkComponentHeader header; 5117b6afd5bSHong Zhang 5127b6afd5bSHong Zhang PetscFunctionBegin; 5137b6afd5bSHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 5147b6afd5bSHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 515e85e6aecSHong Zhang *index = header->index; 5167b6afd5bSHong Zhang PetscFunctionReturn(0); 5177b6afd5bSHong Zhang } 5187b6afd5bSHong Zhang 5195f2c45f1SShri Abhyankar /*@ 520e85e6aecSHong Zhang DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex. 521e85e6aecSHong Zhang 522e85e6aecSHong Zhang Not Collective 523e85e6aecSHong Zhang 524e85e6aecSHong Zhang Input Parameters: 525e85e6aecSHong Zhang + dm - DMNetwork object 526e85e6aecSHong Zhang - p - vertex point 527e85e6aecSHong Zhang 528e85e6aecSHong Zhang Output Paramters: 529e85e6aecSHong Zhang . index - user global numbering for the vertex 530e85e6aecSHong Zhang 531e85e6aecSHong Zhang Level: intermediate 532e85e6aecSHong Zhang 533e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex 534e85e6aecSHong Zhang @*/ 535e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index) 536e85e6aecSHong Zhang { 537e85e6aecSHong Zhang PetscErrorCode ierr; 538e85e6aecSHong Zhang DM_Network *network = (DM_Network*)dm->data; 539e85e6aecSHong Zhang PetscInt offsetp; 540e85e6aecSHong Zhang DMNetworkComponentHeader header; 541e85e6aecSHong Zhang 542e85e6aecSHong Zhang PetscFunctionBegin; 543e85e6aecSHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 544e85e6aecSHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 545e85e6aecSHong Zhang *index = header->index; 546e85e6aecSHong Zhang PetscFunctionReturn(0); 547e85e6aecSHong Zhang } 548e85e6aecSHong Zhang 549c3b11c7cSShri Abhyankar /* 550c3b11c7cSShri Abhyankar DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the 551c3b11c7cSShri Abhyankar component value from the component data array 552c3b11c7cSShri Abhyankar 553c3b11c7cSShri Abhyankar Not Collective 554c3b11c7cSShri Abhyankar 555c3b11c7cSShri Abhyankar Input Parameters: 556c3b11c7cSShri Abhyankar + dm - The DMNetwork object 557c3b11c7cSShri Abhyankar . p - vertex/edge point 558c3b11c7cSShri Abhyankar - compnum - component number 559c3b11c7cSShri Abhyankar 560c3b11c7cSShri Abhyankar Output Parameters: 561c3b11c7cSShri Abhyankar + compkey - the key obtained when registering the component 562c3b11c7cSShri Abhyankar - offset - offset into the component data array associated with the vertex/edge point 563c3b11c7cSShri Abhyankar 564c3b11c7cSShri Abhyankar Notes: 565c3b11c7cSShri Abhyankar Typical usage: 566c3b11c7cSShri Abhyankar 567c3b11c7cSShri Abhyankar DMNetworkGetComponentDataArray(dm, &arr); 568c3b11c7cSShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 569c3b11c7cSShri Abhyankar Loop over vertices or edges 570c3b11c7cSShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 571c3b11c7cSShri Abhyankar Loop over numcomps 572c3b11c7cSShri Abhyankar DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset); 573c3b11c7cSShri Abhyankar compdata = (UserCompDataType)(arr+offset); 574c3b11c7cSShri Abhyankar 575c3b11c7cSShri Abhyankar Level: intermediate 576c3b11c7cSShri Abhyankar 577c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 578c3b11c7cSShri Abhyankar */ 579c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 580c3b11c7cSShri Abhyankar { 581c3b11c7cSShri Abhyankar PetscErrorCode ierr; 582c3b11c7cSShri Abhyankar PetscInt offsetp; 583c3b11c7cSShri Abhyankar DMNetworkComponentHeader header; 584c3b11c7cSShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 585c3b11c7cSShri Abhyankar 586c3b11c7cSShri Abhyankar PetscFunctionBegin; 587c3b11c7cSShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 588c3b11c7cSShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 589c3b11c7cSShri Abhyankar if (compkey) *compkey = header->key[compnum]; 590c3b11c7cSShri Abhyankar if (offset) *offset = offsetp+network->dataheadersize+header->offset[compnum]; 591c3b11c7cSShri Abhyankar PetscFunctionReturn(0); 592c3b11c7cSShri Abhyankar } 593c3b11c7cSShri Abhyankar 594c3b11c7cSShri Abhyankar /*@ 595c3b11c7cSShri Abhyankar DMNetworkGetComponent - Returns the network component and its key 596c3b11c7cSShri Abhyankar 597c3b11c7cSShri Abhyankar Not Collective 598c3b11c7cSShri Abhyankar 599c3b11c7cSShri Abhyankar Input Parameters 600c3b11c7cSShri Abhyankar + dm - DMNetwork object 601c3b11c7cSShri Abhyankar . p - edge or vertex point 602c3b11c7cSShri Abhyankar - compnum - component number 603c3b11c7cSShri Abhyankar 604c3b11c7cSShri Abhyankar Output Parameters: 605c3b11c7cSShri Abhyankar + compkey - the key set for this computing during registration 606c3b11c7cSShri Abhyankar - component - the component data 607c3b11c7cSShri Abhyankar 608c3b11c7cSShri Abhyankar Notes: 609c3b11c7cSShri Abhyankar Typical usage: 610c3b11c7cSShri Abhyankar 611c3b11c7cSShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 612c3b11c7cSShri Abhyankar Loop over vertices or edges 613c3b11c7cSShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 614c3b11c7cSShri Abhyankar Loop over numcomps 615c3b11c7cSShri Abhyankar DMNetworkGetComponent(dm,v,compnum,&key,&component); 616c3b11c7cSShri Abhyankar 617c3b11c7cSShri Abhyankar Level: intermediate 618c3b11c7cSShri Abhyankar 619c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset 620c3b11c7cSShri Abhyankar @*/ 621c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component) 622c3b11c7cSShri Abhyankar { 623c3b11c7cSShri Abhyankar PetscErrorCode ierr; 624c3b11c7cSShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 625c3b11c7cSShri Abhyankar PetscInt offsetd; 626c3b11c7cSShri Abhyankar 627c3b11c7cSShri Abhyankar PetscFunctionBegin; 628c3b11c7cSShri Abhyankar 629c3b11c7cSShri Abhyankar ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr); 630c3b11c7cSShri Abhyankar *component = network->componentdataarray+offsetd; 631c3b11c7cSShri Abhyankar 632c3b11c7cSShri Abhyankar PetscFunctionReturn(0); 633c3b11c7cSShri Abhyankar } 634c3b11c7cSShri Abhyankar 635e85e6aecSHong Zhang /*@ 636325661f6SSatish Balay DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 6375f2c45f1SShri Abhyankar 6385f2c45f1SShri Abhyankar Not Collective 6395f2c45f1SShri Abhyankar 6405f2c45f1SShri Abhyankar Input Parameters: 6415f2c45f1SShri Abhyankar + dm - The DMNetwork object 6425f2c45f1SShri Abhyankar . p - vertex/edge point 6435f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component 6445f2c45f1SShri Abhyankar - compvalue - pointer to the data structure for the component 6455f2c45f1SShri Abhyankar 6465f2c45f1SShri Abhyankar Level: intermediate 6475f2c45f1SShri Abhyankar 6485f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 6495f2c45f1SShri Abhyankar @*/ 6505f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 6515f2c45f1SShri Abhyankar { 6525f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 65343a39a44SBarry Smith DMNetworkComponent *component = &network->component[componentkey]; 6545f2c45f1SShri Abhyankar DMNetworkComponentHeader header = &network->header[p]; 6555f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue = &network->cvalue[p]; 6565f2c45f1SShri Abhyankar PetscErrorCode ierr; 6575f2c45f1SShri Abhyankar 6585f2c45f1SShri Abhyankar PetscFunctionBegin; 659fa58f0a9SHong 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); 660fa58f0a9SHong Zhang 66143a39a44SBarry Smith header->size[header->ndata] = component->size; 66243a39a44SBarry Smith ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 6635f2c45f1SShri Abhyankar header->key[header->ndata] = componentkey; 6645f2c45f1SShri Abhyankar if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 6655f2c45f1SShri Abhyankar 6665f2c45f1SShri Abhyankar cvalue->data[header->ndata] = (void*)compvalue; 6675f2c45f1SShri Abhyankar header->ndata++; 6685f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6695f2c45f1SShri Abhyankar } 6705f2c45f1SShri Abhyankar 6715f2c45f1SShri Abhyankar /*@ 6725f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 6735f2c45f1SShri Abhyankar 6745f2c45f1SShri Abhyankar Not Collective 6755f2c45f1SShri Abhyankar 6765f2c45f1SShri Abhyankar Input Parameters: 6775f2c45f1SShri Abhyankar + dm - The DMNetwork object 6785f2c45f1SShri Abhyankar . p - vertex/edge point 6795f2c45f1SShri Abhyankar 6805f2c45f1SShri Abhyankar Output Parameters: 6815f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 6825f2c45f1SShri Abhyankar 6835f2c45f1SShri Abhyankar Level: intermediate 6845f2c45f1SShri Abhyankar 6855f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 6865f2c45f1SShri Abhyankar @*/ 6875f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 6885f2c45f1SShri Abhyankar { 6895f2c45f1SShri Abhyankar PetscErrorCode ierr; 6905f2c45f1SShri Abhyankar PetscInt offset; 6915f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 6925f2c45f1SShri Abhyankar 6935f2c45f1SShri Abhyankar PetscFunctionBegin; 6945f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 6955f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 6965f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6975f2c45f1SShri Abhyankar } 6985f2c45f1SShri Abhyankar 6995f2c45f1SShri Abhyankar /*@ 7005f2c45f1SShri Abhyankar DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 7015f2c45f1SShri Abhyankar 7025f2c45f1SShri Abhyankar Not Collective 7035f2c45f1SShri Abhyankar 7045f2c45f1SShri Abhyankar Input Parameters: 7055f2c45f1SShri Abhyankar + dm - The DMNetwork object 7065f2c45f1SShri Abhyankar - p - the edge/vertex point 7075f2c45f1SShri Abhyankar 7085f2c45f1SShri Abhyankar Output Parameters: 7095f2c45f1SShri Abhyankar . offset - the offset 7105f2c45f1SShri Abhyankar 7115f2c45f1SShri Abhyankar Level: intermediate 7125f2c45f1SShri Abhyankar 7135f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 7145f2c45f1SShri Abhyankar @*/ 7155f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 7165f2c45f1SShri Abhyankar { 7175f2c45f1SShri Abhyankar PetscErrorCode ierr; 7185f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7195f2c45f1SShri Abhyankar 7205f2c45f1SShri Abhyankar PetscFunctionBegin; 7215f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr); 7225f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7235f2c45f1SShri Abhyankar } 7245f2c45f1SShri Abhyankar 7255f2c45f1SShri Abhyankar /*@ 7265f2c45f1SShri Abhyankar DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 7275f2c45f1SShri Abhyankar 7285f2c45f1SShri Abhyankar Not Collective 7295f2c45f1SShri Abhyankar 7305f2c45f1SShri Abhyankar Input Parameters: 7315f2c45f1SShri Abhyankar + dm - The DMNetwork object 7325f2c45f1SShri Abhyankar - p - the edge/vertex point 7335f2c45f1SShri Abhyankar 7345f2c45f1SShri Abhyankar Output Parameters: 7355f2c45f1SShri Abhyankar . offsetg - the offset 7365f2c45f1SShri Abhyankar 7375f2c45f1SShri Abhyankar Level: intermediate 7385f2c45f1SShri Abhyankar 7395f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 7405f2c45f1SShri Abhyankar @*/ 7415f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 7425f2c45f1SShri Abhyankar { 7435f2c45f1SShri Abhyankar PetscErrorCode ierr; 7445f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7455f2c45f1SShri Abhyankar 7465f2c45f1SShri Abhyankar PetscFunctionBegin; 7475f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr); 7486fefedf4SHong Zhang if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */ 7495f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7505f2c45f1SShri Abhyankar } 7515f2c45f1SShri Abhyankar 75224121865SAdrian Maldonado /*@ 75324121865SAdrian Maldonado DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector. 75424121865SAdrian Maldonado 75524121865SAdrian Maldonado Not Collective 75624121865SAdrian Maldonado 75724121865SAdrian Maldonado Input Parameters: 75824121865SAdrian Maldonado + dm - The DMNetwork object 75924121865SAdrian Maldonado - p - the edge point 76024121865SAdrian Maldonado 76124121865SAdrian Maldonado Output Parameters: 76224121865SAdrian Maldonado . offset - the offset 76324121865SAdrian Maldonado 76424121865SAdrian Maldonado Level: intermediate 76524121865SAdrian Maldonado 76624121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 76724121865SAdrian Maldonado @*/ 76824121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset) 76924121865SAdrian Maldonado { 77024121865SAdrian Maldonado PetscErrorCode ierr; 77124121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 77224121865SAdrian Maldonado 77324121865SAdrian Maldonado PetscFunctionBegin; 77424121865SAdrian Maldonado 77524121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr); 77624121865SAdrian Maldonado PetscFunctionReturn(0); 77724121865SAdrian Maldonado } 77824121865SAdrian Maldonado 77924121865SAdrian Maldonado /*@ 78024121865SAdrian Maldonado DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector. 78124121865SAdrian Maldonado 78224121865SAdrian Maldonado Not Collective 78324121865SAdrian Maldonado 78424121865SAdrian Maldonado Input Parameters: 78524121865SAdrian Maldonado + dm - The DMNetwork object 78624121865SAdrian Maldonado - p - the vertex point 78724121865SAdrian Maldonado 78824121865SAdrian Maldonado Output Parameters: 78924121865SAdrian Maldonado . offset - the offset 79024121865SAdrian Maldonado 79124121865SAdrian Maldonado Level: intermediate 79224121865SAdrian Maldonado 79324121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 79424121865SAdrian Maldonado @*/ 79524121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset) 79624121865SAdrian Maldonado { 79724121865SAdrian Maldonado PetscErrorCode ierr; 79824121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 79924121865SAdrian Maldonado 80024121865SAdrian Maldonado PetscFunctionBegin; 80124121865SAdrian Maldonado 80224121865SAdrian Maldonado p -= network->vStart; 80324121865SAdrian Maldonado 80424121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr); 80524121865SAdrian Maldonado PetscFunctionReturn(0); 80624121865SAdrian Maldonado } 8075f2c45f1SShri Abhyankar /*@ 8085f2c45f1SShri Abhyankar DMNetworkAddNumVariables - Add number of variables associated with a given point. 8095f2c45f1SShri Abhyankar 8105f2c45f1SShri Abhyankar Not Collective 8115f2c45f1SShri Abhyankar 8125f2c45f1SShri Abhyankar Input Parameters: 8135f2c45f1SShri Abhyankar + dm - The DMNetworkObject 8145f2c45f1SShri Abhyankar . p - the vertex/edge point 8155f2c45f1SShri Abhyankar - nvar - number of additional variables 8165f2c45f1SShri Abhyankar 8175f2c45f1SShri Abhyankar Level: intermediate 8185f2c45f1SShri Abhyankar 8195f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables 8205f2c45f1SShri Abhyankar @*/ 8215f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 8225f2c45f1SShri Abhyankar { 8235f2c45f1SShri Abhyankar PetscErrorCode ierr; 8245f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8255f2c45f1SShri Abhyankar 8265f2c45f1SShri Abhyankar PetscFunctionBegin; 8275f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 8285f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8295f2c45f1SShri Abhyankar } 8305f2c45f1SShri Abhyankar 83127f51fceSHong Zhang /*@ 83227f51fceSHong Zhang DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point. 83327f51fceSHong Zhang 83427f51fceSHong Zhang Not Collective 83527f51fceSHong Zhang 83627f51fceSHong Zhang Input Parameters: 83727f51fceSHong Zhang + dm - The DMNetworkObject 83827f51fceSHong Zhang - p - the vertex/edge point 83927f51fceSHong Zhang 84027f51fceSHong Zhang Output Parameters: 84127f51fceSHong Zhang . nvar - number of variables 84227f51fceSHong Zhang 84327f51fceSHong Zhang Level: intermediate 84427f51fceSHong Zhang 84527f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables 84627f51fceSHong Zhang @*/ 84727f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar) 84827f51fceSHong Zhang { 84927f51fceSHong Zhang PetscErrorCode ierr; 85027f51fceSHong Zhang DM_Network *network = (DM_Network*)dm->data; 85127f51fceSHong Zhang 85227f51fceSHong Zhang PetscFunctionBegin; 85327f51fceSHong Zhang ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 85427f51fceSHong Zhang PetscFunctionReturn(0); 85527f51fceSHong Zhang } 85627f51fceSHong Zhang 8575f2c45f1SShri Abhyankar /*@ 8585f2c45f1SShri Abhyankar DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 8595f2c45f1SShri Abhyankar 8605f2c45f1SShri Abhyankar Not Collective 8615f2c45f1SShri Abhyankar 8625f2c45f1SShri Abhyankar Input Parameters: 8635f2c45f1SShri Abhyankar + dm - The DMNetworkObject 8645f2c45f1SShri Abhyankar . p - the vertex/edge point 8655f2c45f1SShri Abhyankar - nvar - number of variables 8665f2c45f1SShri Abhyankar 8675f2c45f1SShri Abhyankar Level: intermediate 8685f2c45f1SShri Abhyankar 8695f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables 8705f2c45f1SShri Abhyankar @*/ 8715f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 8725f2c45f1SShri Abhyankar { 8735f2c45f1SShri Abhyankar PetscErrorCode ierr; 8745f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8755f2c45f1SShri Abhyankar 8765f2c45f1SShri Abhyankar PetscFunctionBegin; 8775f2c45f1SShri Abhyankar ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 8785f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8795f2c45f1SShri Abhyankar } 8805f2c45f1SShri Abhyankar 8815f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This 8825f2c45f1SShri Abhyankar function is called during DMSetUp() */ 8835f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 8845f2c45f1SShri Abhyankar { 8855f2c45f1SShri Abhyankar PetscErrorCode ierr; 8865f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8875f2c45f1SShri Abhyankar PetscInt arr_size; 8885f2c45f1SShri Abhyankar PetscInt p,offset,offsetp; 8895f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 8905f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 8915f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 8925f2c45f1SShri Abhyankar PetscInt ncomp, i; 8935f2c45f1SShri Abhyankar 8945f2c45f1SShri Abhyankar PetscFunctionBegin; 8955f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 8965f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 89775b160a0SShri Abhyankar ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 8985f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 8995f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 9005f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 9015f2c45f1SShri Abhyankar /* Copy header */ 9025f2c45f1SShri Abhyankar header = &network->header[p]; 903302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 9045f2c45f1SShri Abhyankar /* Copy data */ 9055f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 9065f2c45f1SShri Abhyankar ncomp = header->ndata; 9075f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 9085f2c45f1SShri Abhyankar offset = offsetp + network->dataheadersize + header->offset[i]; 909302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 9105f2c45f1SShri Abhyankar } 9115f2c45f1SShri Abhyankar } 9125f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9135f2c45f1SShri Abhyankar } 9145f2c45f1SShri Abhyankar 9155f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 9165f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm) 9175f2c45f1SShri Abhyankar { 9185f2c45f1SShri Abhyankar PetscErrorCode ierr; 9195f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9205f2c45f1SShri Abhyankar 9215f2c45f1SShri Abhyankar PetscFunctionBegin; 9225f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 9235f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9245f2c45f1SShri Abhyankar } 9255f2c45f1SShri Abhyankar 9265f2c45f1SShri Abhyankar /*@C 9275f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray - Returns the component data array 9285f2c45f1SShri Abhyankar 9295f2c45f1SShri Abhyankar Not Collective 9305f2c45f1SShri Abhyankar 9315f2c45f1SShri Abhyankar Input Parameters: 9325f2c45f1SShri Abhyankar . dm - The DMNetwork Object 9335f2c45f1SShri Abhyankar 9345f2c45f1SShri Abhyankar Output Parameters: 9355f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components 9365f2c45f1SShri Abhyankar 9375f2c45f1SShri Abhyankar Level: intermediate 9385f2c45f1SShri Abhyankar 939a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents 9405f2c45f1SShri Abhyankar @*/ 9415f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 9425f2c45f1SShri Abhyankar { 9435f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9445f2c45f1SShri Abhyankar 9455f2c45f1SShri Abhyankar PetscFunctionBegin; 9465f2c45f1SShri Abhyankar *componentdataarray = network->componentdataarray; 9475f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9485f2c45f1SShri Abhyankar } 9495f2c45f1SShri Abhyankar 95024121865SAdrian Maldonado /* Get a subsection from a range of points */ 95124121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection) 95224121865SAdrian Maldonado { 95324121865SAdrian Maldonado PetscErrorCode ierr; 95424121865SAdrian Maldonado PetscInt i, nvar; 95524121865SAdrian Maldonado 95624121865SAdrian Maldonado PetscFunctionBegin; 95724121865SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr); 95824121865SAdrian Maldonado ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr); 95924121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 96024121865SAdrian Maldonado ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr); 96124121865SAdrian Maldonado ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr); 96224121865SAdrian Maldonado } 96324121865SAdrian Maldonado 96424121865SAdrian Maldonado ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr); 96524121865SAdrian Maldonado PetscFunctionReturn(0); 96624121865SAdrian Maldonado } 96724121865SAdrian Maldonado 96824121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */ 96924121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 97024121865SAdrian Maldonado { 97124121865SAdrian Maldonado PetscErrorCode ierr; 97224121865SAdrian Maldonado PetscInt i, *subpoints; 97324121865SAdrian Maldonado 97424121865SAdrian Maldonado PetscFunctionBegin; 97524121865SAdrian Maldonado /* Create index sets to map from "points" to "subpoints" */ 97624121865SAdrian Maldonado ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr); 97724121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 97824121865SAdrian Maldonado subpoints[i - pstart] = i; 97924121865SAdrian Maldonado } 980459726d8SSatish Balay ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr); 98124121865SAdrian Maldonado ierr = PetscFree(subpoints);CHKERRQ(ierr); 98224121865SAdrian Maldonado PetscFunctionReturn(0); 98324121865SAdrian Maldonado } 98424121865SAdrian Maldonado 98524121865SAdrian Maldonado /*@ 98624121865SAdrian Maldonado DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute. 98724121865SAdrian Maldonado 98824121865SAdrian Maldonado Collective 98924121865SAdrian Maldonado 99024121865SAdrian Maldonado Input Parameters: 99124121865SAdrian Maldonado . dm - The DMNetworkObject 99224121865SAdrian Maldonado 99324121865SAdrian Maldonado Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 99424121865SAdrian Maldonado 99524121865SAdrian Maldonado points = [0 1 2 3 4 5 6] 99624121865SAdrian Maldonado 99724121865SAdrian Maldonado where edges = [0, 3] and vertices = [4, 6]. The new orderings will be specific to the subset (i.e vertices = [0, 2]). 99824121865SAdrian Maldonado 99924121865SAdrian Maldonado With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 100024121865SAdrian Maldonado 100124121865SAdrian Maldonado Level: intermediate 100224121865SAdrian Maldonado 100324121865SAdrian Maldonado @*/ 100424121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 100524121865SAdrian Maldonado { 100624121865SAdrian Maldonado PetscErrorCode ierr; 100724121865SAdrian Maldonado MPI_Comm comm; 10089852e123SBarry Smith PetscMPIInt rank, size; 100924121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 101024121865SAdrian Maldonado 1011eab1376dSHong Zhang PetscFunctionBegin; 101224121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 101324121865SAdrian Maldonado ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 10149852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 101524121865SAdrian Maldonado 101624121865SAdrian Maldonado /* Create maps for vertices and edges */ 101724121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 101824121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr); 101924121865SAdrian Maldonado 102024121865SAdrian Maldonado /* Create local sub-sections */ 102124121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr); 102224121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr); 102324121865SAdrian Maldonado 10249852e123SBarry Smith if (size > 1) { 102524121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 102624121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr); 102724121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr); 102824121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr); 102924121865SAdrian Maldonado } else { 103024121865SAdrian Maldonado /* create structures for vertex */ 103124121865SAdrian Maldonado ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr); 103224121865SAdrian Maldonado /* create structures for edge */ 103324121865SAdrian Maldonado ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr); 103424121865SAdrian Maldonado } 103524121865SAdrian Maldonado 103624121865SAdrian Maldonado 103724121865SAdrian Maldonado /* Add viewers */ 103824121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr); 103924121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr); 104024121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr); 104124121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr); 104224121865SAdrian Maldonado 104324121865SAdrian Maldonado PetscFunctionReturn(0); 104424121865SAdrian Maldonado } 10457b6afd5bSHong Zhang 10465f2c45f1SShri Abhyankar /*@ 10475f2c45f1SShri Abhyankar DMNetworkDistribute - Distributes the network and moves associated component data. 10485f2c45f1SShri Abhyankar 10495f2c45f1SShri Abhyankar Collective 10505f2c45f1SShri Abhyankar 10515f2c45f1SShri Abhyankar Input Parameter: 1052d3464fd4SAdrian Maldonado + DM - the DMNetwork object 10535f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default 10545f2c45f1SShri Abhyankar 10555f2c45f1SShri Abhyankar Notes: 10568b171c8eSHong Zhang Distributes the network with <overlap>-overlapping partitioning of the edges. 10575f2c45f1SShri Abhyankar 10585f2c45f1SShri Abhyankar Level: intermediate 10595f2c45f1SShri Abhyankar 10605f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 10615f2c45f1SShri Abhyankar @*/ 1062d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap) 10635f2c45f1SShri Abhyankar { 1064d3464fd4SAdrian Maldonado MPI_Comm comm; 10655f2c45f1SShri Abhyankar PetscErrorCode ierr; 1066d3464fd4SAdrian Maldonado PetscMPIInt size; 1067d3464fd4SAdrian Maldonado DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data); 1068d3464fd4SAdrian Maldonado DM_Network *newDMnetwork; 10695f2c45f1SShri Abhyankar PetscSF pointsf; 10705f2c45f1SShri Abhyankar DM newDM; 107151ac5effSHong Zhang PetscPartitioner part; 1072b9c6e19dSShri Abhyankar PetscInt j,e,v,offset; 1073b9c6e19dSShri Abhyankar DMNetworkComponentHeader header; 10745f2c45f1SShri Abhyankar 10755f2c45f1SShri Abhyankar PetscFunctionBegin; 1076d3464fd4SAdrian Maldonado 1077d3464fd4SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr); 1078d3464fd4SAdrian Maldonado ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1079d3464fd4SAdrian Maldonado if (size == 1) PetscFunctionReturn(0); 1080d3464fd4SAdrian Maldonado 1081d3464fd4SAdrian Maldonado ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr); 10825f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 10835f2c45f1SShri Abhyankar newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 108451ac5effSHong Zhang 108551ac5effSHong Zhang /* Enable runtime options for petscpartitioner */ 108651ac5effSHong Zhang ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr); 108751ac5effSHong Zhang ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 108851ac5effSHong Zhang 10895f2c45f1SShri Abhyankar /* Distribute plex dm and dof section */ 109080cf41d5SMatthew G. Knepley ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 109151ac5effSHong Zhang 10925f2c45f1SShri Abhyankar /* Distribute dof section */ 1093d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr); 10945f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 1095d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr); 109651ac5effSHong Zhang 10975f2c45f1SShri Abhyankar /* Distribute data and associated section */ 109831da1fc8SHong Zhang ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 109924121865SAdrian Maldonado 11005f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 11015f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 11025f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 11035f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 11046fefedf4SHong Zhang newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart; 11056fefedf4SHong Zhang newDMnetwork->NVertices = oldDMnetwork->NVertices; 11065f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 110724121865SAdrian Maldonado 11085f2c45f1SShri Abhyankar /* Set Dof section as the default section for dm */ 11095f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 11105f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 11115f2c45f1SShri Abhyankar 1112b9c6e19dSShri Abhyankar /* Set up subnetwork info in the newDM */ 1113b9c6e19dSShri Abhyankar newDMnetwork->nsubnet = oldDMnetwork->nsubnet; 1114b9c6e19dSShri Abhyankar ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr); 1115b9c6e19dSShri Abhyankar /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already 1116b9c6e19dSShri Abhyankar calculated in DMNetworkLayoutSetUp() 1117b9c6e19dSShri Abhyankar */ 1118b9c6e19dSShri Abhyankar for(j=0; j < newDMnetwork->nsubnet; j++) { 1119b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx; 1120b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge; 1121b9c6e19dSShri Abhyankar } 1122b9c6e19dSShri Abhyankar 1123b9c6e19dSShri Abhyankar for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) { 1124b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1125b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1126b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].nedge++; 1127b9c6e19dSShri Abhyankar } 1128b9c6e19dSShri Abhyankar 1129b9c6e19dSShri Abhyankar for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) { 1130b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1131b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1132b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].nvtx++; 1133b9c6e19dSShri Abhyankar } 1134b9c6e19dSShri Abhyankar 1135b9c6e19dSShri Abhyankar /* Now create the vertices and edge arrays for the subnetworks */ 1136b9c6e19dSShri Abhyankar for(j=0; j < newDMnetwork->nsubnet; j++) { 1137b9c6e19dSShri Abhyankar ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr); 1138b9c6e19dSShri Abhyankar ierr = PetscCalloc1(newDMnetwork->subnet[j].nvtx,&newDMnetwork->subnet[j].vertices);CHKERRQ(ierr); 1139b9c6e19dSShri Abhyankar /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 1140b9c6e19dSShri Abhyankar These get updated when the vertices and edges are added. */ 1141b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0; 1142b9c6e19dSShri Abhyankar } 1143b9c6e19dSShri Abhyankar 1144b9c6e19dSShri Abhyankar /* Set the vertices and edges in each subnetwork */ 1145b9c6e19dSShri Abhyankar for(e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) { 1146b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1147b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1148b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e; 1149b9c6e19dSShri Abhyankar } 1150b9c6e19dSShri Abhyankar 1151b9c6e19dSShri Abhyankar for(v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) { 1152b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1153b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1154b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v; 1155b9c6e19dSShri Abhyankar } 1156b9c6e19dSShri Abhyankar 115724121865SAdrian Maldonado /* Destroy point SF */ 115824121865SAdrian Maldonado ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 115924121865SAdrian Maldonado 1160d3464fd4SAdrian Maldonado ierr = DMDestroy(dm);CHKERRQ(ierr); 1161d3464fd4SAdrian Maldonado *dm = newDM; 11625f2c45f1SShri Abhyankar PetscFunctionReturn(0); 11635f2c45f1SShri Abhyankar } 11645f2c45f1SShri Abhyankar 116524121865SAdrian Maldonado /*@C 116624121865SAdrian Maldonado PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering. 116724121865SAdrian Maldonado 116824121865SAdrian Maldonado Input Parameters: 116924121865SAdrian Maldonado + masterSF - the original SF structure 117024121865SAdrian Maldonado - map - a ISLocalToGlobal mapping that contains the subset of points 117124121865SAdrian Maldonado 117224121865SAdrian Maldonado Output Parameters: 117324121865SAdrian Maldonado . subSF - a subset of the masterSF for the desired subset. 117424121865SAdrian Maldonado */ 117524121865SAdrian Maldonado 117624121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) { 117724121865SAdrian Maldonado 117824121865SAdrian Maldonado PetscErrorCode ierr; 117924121865SAdrian Maldonado PetscInt nroots, nleaves, *ilocal_sub; 118024121865SAdrian Maldonado PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 118124121865SAdrian Maldonado PetscInt *local_points, *remote_points; 118224121865SAdrian Maldonado PetscSFNode *iremote_sub; 118324121865SAdrian Maldonado const PetscInt *ilocal; 118424121865SAdrian Maldonado const PetscSFNode *iremote; 118524121865SAdrian Maldonado 118624121865SAdrian Maldonado PetscFunctionBegin; 118724121865SAdrian Maldonado ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 118824121865SAdrian Maldonado 118924121865SAdrian Maldonado /* Look for leaves that pertain to the subset of points. Get the local ordering */ 119024121865SAdrian Maldonado ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr); 119124121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr); 119224121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 119324121865SAdrian Maldonado if (ilocal_map[i] != -1) nleaves_sub += 1; 119424121865SAdrian Maldonado } 119524121865SAdrian Maldonado /* Re-number ilocal with subset numbering. Need information from roots */ 119624121865SAdrian Maldonado ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr); 119724121865SAdrian Maldonado for (i = 0; i < nroots; i++) local_points[i] = i; 119824121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr); 119924121865SAdrian Maldonado ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 120024121865SAdrian Maldonado ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 120124121865SAdrian Maldonado /* Fill up graph using local (that is, local to the subset) numbering. */ 12024b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr); 12034b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr); 120424121865SAdrian Maldonado nleaves_sub = 0; 120524121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 120624121865SAdrian Maldonado if (ilocal_map[i] != -1) { 120724121865SAdrian Maldonado ilocal_sub[nleaves_sub] = ilocal_map[i]; 12084b70a8deSAdrian Maldonado iremote_sub[nleaves_sub].rank = iremote[i].rank; 120924121865SAdrian Maldonado iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 121024121865SAdrian Maldonado nleaves_sub += 1; 121124121865SAdrian Maldonado } 121224121865SAdrian Maldonado } 121324121865SAdrian Maldonado ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr); 121424121865SAdrian Maldonado ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr); 121524121865SAdrian Maldonado 121624121865SAdrian Maldonado /* Create new subSF */ 121724121865SAdrian Maldonado ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr); 121824121865SAdrian Maldonado ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr); 12194b70a8deSAdrian Maldonado ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr); 122024121865SAdrian Maldonado ierr = PetscFree(ilocal_map);CHKERRQ(ierr); 12214b70a8deSAdrian Maldonado ierr = PetscFree(iremote_sub);CHKERRQ(ierr); 122224121865SAdrian Maldonado PetscFunctionReturn(0); 122324121865SAdrian Maldonado } 122424121865SAdrian Maldonado 12255f2c45f1SShri Abhyankar /*@C 12265f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 12275f2c45f1SShri Abhyankar 12285f2c45f1SShri Abhyankar Not Collective 12295f2c45f1SShri Abhyankar 12305f2c45f1SShri Abhyankar Input Parameters: 12315f2c45f1SShri Abhyankar + dm - The DMNetwork object 12325f2c45f1SShri Abhyankar - p - the vertex point 12335f2c45f1SShri Abhyankar 12345f2c45f1SShri Abhyankar Output Paramters: 12355f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 12365f2c45f1SShri Abhyankar - edges - List of edge points 12375f2c45f1SShri Abhyankar 12385f2c45f1SShri Abhyankar Level: intermediate 12395f2c45f1SShri Abhyankar 12405f2c45f1SShri Abhyankar Fortran Notes: 12415f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 12425f2c45f1SShri Abhyankar include petsc.h90 in your code. 12435f2c45f1SShri Abhyankar 1244d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices 12455f2c45f1SShri Abhyankar @*/ 12465f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 12475f2c45f1SShri Abhyankar { 12485f2c45f1SShri Abhyankar PetscErrorCode ierr; 12495f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 12505f2c45f1SShri Abhyankar 12515f2c45f1SShri Abhyankar PetscFunctionBegin; 12525f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 12535f2c45f1SShri Abhyankar ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 12545f2c45f1SShri Abhyankar PetscFunctionReturn(0); 12555f2c45f1SShri Abhyankar } 12565f2c45f1SShri Abhyankar 12575f2c45f1SShri Abhyankar /*@C 1258d842c372SHong Zhang DMNetworkGetConnectedVertices - Return the connected vertices for this edge point 12595f2c45f1SShri Abhyankar 12605f2c45f1SShri Abhyankar Not Collective 12615f2c45f1SShri Abhyankar 12625f2c45f1SShri Abhyankar Input Parameters: 12635f2c45f1SShri Abhyankar + dm - The DMNetwork object 12645f2c45f1SShri Abhyankar - p - the edge point 12655f2c45f1SShri Abhyankar 12665f2c45f1SShri Abhyankar Output Paramters: 12675f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 12685f2c45f1SShri Abhyankar 12695f2c45f1SShri Abhyankar Level: intermediate 12705f2c45f1SShri Abhyankar 12715f2c45f1SShri Abhyankar Fortran Notes: 12725f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 12735f2c45f1SShri Abhyankar include petsc.h90 in your code. 12745f2c45f1SShri Abhyankar 12755f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 12765f2c45f1SShri Abhyankar @*/ 1277d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[]) 12785f2c45f1SShri Abhyankar { 12795f2c45f1SShri Abhyankar PetscErrorCode ierr; 12805f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 12815f2c45f1SShri Abhyankar 12825f2c45f1SShri Abhyankar PetscFunctionBegin; 12835f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 12845f2c45f1SShri Abhyankar PetscFunctionReturn(0); 12855f2c45f1SShri Abhyankar } 12865f2c45f1SShri Abhyankar 12875f2c45f1SShri Abhyankar /*@ 12885f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 12895f2c45f1SShri Abhyankar 12905f2c45f1SShri Abhyankar Not Collective 12915f2c45f1SShri Abhyankar 12925f2c45f1SShri Abhyankar Input Parameters: 12935f2c45f1SShri Abhyankar + dm - The DMNetwork object 12945f2c45f1SShri Abhyankar . p - the vertex point 12955f2c45f1SShri Abhyankar 12965f2c45f1SShri Abhyankar Output Parameter: 12975f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 12985f2c45f1SShri Abhyankar 12995f2c45f1SShri Abhyankar Level: intermediate 13005f2c45f1SShri Abhyankar 1301d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange 13025f2c45f1SShri Abhyankar @*/ 13035f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 13045f2c45f1SShri Abhyankar { 13055f2c45f1SShri Abhyankar PetscErrorCode ierr; 13065f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 13075f2c45f1SShri Abhyankar PetscInt offsetg; 13085f2c45f1SShri Abhyankar PetscSection sectiong; 13095f2c45f1SShri Abhyankar 13105f2c45f1SShri Abhyankar PetscFunctionBegin; 13115f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 13125f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,§iong);CHKERRQ(ierr); 13135f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 13145f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 13155f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13165f2c45f1SShri Abhyankar } 13175f2c45f1SShri Abhyankar 13185f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 13195f2c45f1SShri Abhyankar { 13205f2c45f1SShri Abhyankar PetscErrorCode ierr; 13215f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 13225f2c45f1SShri Abhyankar 13235f2c45f1SShri Abhyankar PetscFunctionBegin; 13245f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 13255f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 13265f2c45f1SShri Abhyankar 13275f2c45f1SShri Abhyankar ierr = DMSetDefaultSection(network->plex,network->DofSection);CHKERRQ(ierr); 13285f2c45f1SShri Abhyankar ierr = DMGetDefaultGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 13295f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13305f2c45f1SShri Abhyankar } 13315f2c45f1SShri Abhyankar 13321ad426b7SHong Zhang /*@ 133317df6e9eSHong Zhang DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 13341ad426b7SHong Zhang -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 13351ad426b7SHong Zhang 13361ad426b7SHong Zhang Collective 13371ad426b7SHong Zhang 13381ad426b7SHong Zhang Input Parameters: 133983b2e829SHong Zhang + dm - The DMNetwork object 134083b2e829SHong Zhang . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 134183b2e829SHong Zhang - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 13421ad426b7SHong Zhang 13431ad426b7SHong Zhang Level: intermediate 13441ad426b7SHong Zhang 13451ad426b7SHong Zhang @*/ 134683b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 13471ad426b7SHong Zhang { 13481ad426b7SHong Zhang DM_Network *network=(DM_Network*)dm->data; 13498675203cSHong Zhang PetscErrorCode ierr; 13501ad426b7SHong Zhang 13511ad426b7SHong Zhang PetscFunctionBegin; 135283b2e829SHong Zhang network->userEdgeJacobian = eflg; 135383b2e829SHong Zhang network->userVertexJacobian = vflg; 13548675203cSHong Zhang 13558675203cSHong Zhang if (eflg && !network->Je) { 13568675203cSHong Zhang ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 13578675203cSHong Zhang } 13588675203cSHong Zhang 13598675203cSHong Zhang if (vflg && !network->Jv) { 13608675203cSHong Zhang PetscInt i,*vptr,nedges,vStart=network->vStart; 13618675203cSHong Zhang PetscInt nVertices = network->nVertices,nedges_total; 13628675203cSHong Zhang const PetscInt *edges; 13638675203cSHong Zhang 13648675203cSHong Zhang /* count nvertex_total */ 13658675203cSHong Zhang nedges_total = 0; 13668675203cSHong Zhang ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr); 13678675203cSHong Zhang 13688675203cSHong Zhang vptr[0] = 0; 13698675203cSHong Zhang for (i=0; i<nVertices; i++) { 13708675203cSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 13718675203cSHong Zhang nedges_total += nedges; 13728675203cSHong Zhang vptr[i+1] = vptr[i] + 2*nedges + 1; 13738675203cSHong Zhang } 13748675203cSHong Zhang 13758675203cSHong Zhang ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr); 13768675203cSHong Zhang network->Jvptr = vptr; 13778675203cSHong Zhang } 13781ad426b7SHong Zhang PetscFunctionReturn(0); 13791ad426b7SHong Zhang } 13801ad426b7SHong Zhang 13811ad426b7SHong Zhang /*@ 138283b2e829SHong Zhang DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 138383b2e829SHong Zhang 138483b2e829SHong Zhang Not Collective 138583b2e829SHong Zhang 138683b2e829SHong Zhang Input Parameters: 138783b2e829SHong Zhang + dm - The DMNetwork object 138883b2e829SHong Zhang . p - the edge point 13893e97b6e8SHong Zhang - J - array (size = 3) of Jacobian submatrices for this edge point: 13903e97b6e8SHong Zhang J[0]: this edge 1391d842c372SHong Zhang J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices() 139283b2e829SHong Zhang 139383b2e829SHong Zhang Level: intermediate 139483b2e829SHong Zhang 139583b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix 139683b2e829SHong Zhang @*/ 139783b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 139883b2e829SHong Zhang { 139983b2e829SHong Zhang DM_Network *network=(DM_Network*)dm->data; 140083b2e829SHong Zhang 140183b2e829SHong Zhang PetscFunctionBegin; 14028675203cSHong Zhang if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 14038675203cSHong Zhang 14048675203cSHong Zhang if (J) { 1405883e35e8SHong Zhang network->Je[3*p] = J[0]; 1406883e35e8SHong Zhang network->Je[3*p+1] = J[1]; 1407883e35e8SHong Zhang network->Je[3*p+2] = J[2]; 14088675203cSHong Zhang } 140983b2e829SHong Zhang PetscFunctionReturn(0); 141083b2e829SHong Zhang } 141183b2e829SHong Zhang 141283b2e829SHong Zhang /*@ 141376ddfea5SHong Zhang DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 14141ad426b7SHong Zhang 14151ad426b7SHong Zhang Not Collective 14161ad426b7SHong Zhang 14171ad426b7SHong Zhang Input Parameters: 14181ad426b7SHong Zhang + dm - The DMNetwork object 14191ad426b7SHong Zhang . p - the vertex point 14203e97b6e8SHong Zhang - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 14213e97b6e8SHong Zhang J[0]: this vertex 14223e97b6e8SHong Zhang J[1+2*i]: i-th supporting edge 14233e97b6e8SHong Zhang J[1+2*i+1]: i-th connected vertex 14241ad426b7SHong Zhang 14251ad426b7SHong Zhang Level: intermediate 14261ad426b7SHong Zhang 142783b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix 14281ad426b7SHong Zhang @*/ 1429883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 14305f2c45f1SShri Abhyankar { 14315f2c45f1SShri Abhyankar PetscErrorCode ierr; 14325f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 14338675203cSHong Zhang PetscInt i,*vptr,nedges,vStart=network->vStart; 1434883e35e8SHong Zhang const PetscInt *edges; 14355f2c45f1SShri Abhyankar 14365f2c45f1SShri Abhyankar PetscFunctionBegin; 14378675203cSHong Zhang if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 1438883e35e8SHong Zhang 14398675203cSHong Zhang if (J) { 1440883e35e8SHong Zhang vptr = network->Jvptr; 14413e97b6e8SHong Zhang network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 14423e97b6e8SHong Zhang 14433e97b6e8SHong Zhang /* Set Jacobian for each supporting edge and connected vertex */ 1444883e35e8SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 1445883e35e8SHong Zhang for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 14468675203cSHong Zhang } 1447883e35e8SHong Zhang PetscFunctionReturn(0); 1448883e35e8SHong Zhang } 1449883e35e8SHong Zhang 1450e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 14515cf7da58SHong Zhang { 14525cf7da58SHong Zhang PetscErrorCode ierr; 14535cf7da58SHong Zhang PetscInt j; 14545cf7da58SHong Zhang PetscScalar val=(PetscScalar)ncols; 14555cf7da58SHong Zhang 14565cf7da58SHong Zhang PetscFunctionBegin; 14575cf7da58SHong Zhang if (!ghost) { 14585cf7da58SHong Zhang for (j=0; j<nrows; j++) { 14595cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 14605cf7da58SHong Zhang } 14615cf7da58SHong Zhang } else { 14625cf7da58SHong Zhang for (j=0; j<nrows; j++) { 14635cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 14645cf7da58SHong Zhang } 14655cf7da58SHong Zhang } 14665cf7da58SHong Zhang PetscFunctionReturn(0); 14675cf7da58SHong Zhang } 14685cf7da58SHong Zhang 1469e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 14705cf7da58SHong Zhang { 14715cf7da58SHong Zhang PetscErrorCode ierr; 14725cf7da58SHong Zhang PetscInt j,ncols_u; 14735cf7da58SHong Zhang PetscScalar val; 14745cf7da58SHong Zhang 14755cf7da58SHong Zhang PetscFunctionBegin; 14765cf7da58SHong Zhang if (!ghost) { 14775cf7da58SHong Zhang for (j=0; j<nrows; j++) { 14785cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 14795cf7da58SHong Zhang val = (PetscScalar)ncols_u; 14805cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 14815cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 14825cf7da58SHong Zhang } 14835cf7da58SHong Zhang } else { 14845cf7da58SHong Zhang for (j=0; j<nrows; j++) { 14855cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 14865cf7da58SHong Zhang val = (PetscScalar)ncols_u; 14875cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 14885cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 14895cf7da58SHong Zhang } 14905cf7da58SHong Zhang } 14915cf7da58SHong Zhang PetscFunctionReturn(0); 14925cf7da58SHong Zhang } 14935cf7da58SHong Zhang 1494e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 14955cf7da58SHong Zhang { 14965cf7da58SHong Zhang PetscErrorCode ierr; 14975cf7da58SHong Zhang 14985cf7da58SHong Zhang PetscFunctionBegin; 14995cf7da58SHong Zhang if (Ju) { 15005cf7da58SHong Zhang ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 15015cf7da58SHong Zhang } else { 15025cf7da58SHong Zhang ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 15035cf7da58SHong Zhang } 15045cf7da58SHong Zhang PetscFunctionReturn(0); 15055cf7da58SHong Zhang } 15065cf7da58SHong Zhang 1507e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1508883e35e8SHong Zhang { 1509883e35e8SHong Zhang PetscErrorCode ierr; 1510883e35e8SHong Zhang PetscInt j,*cols; 1511883e35e8SHong Zhang PetscScalar *zeros; 1512883e35e8SHong Zhang 1513883e35e8SHong Zhang PetscFunctionBegin; 1514883e35e8SHong Zhang ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 1515883e35e8SHong Zhang for (j=0; j<ncols; j++) cols[j] = j+ cstart; 1516883e35e8SHong Zhang ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 1517883e35e8SHong Zhang ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 15181ad426b7SHong Zhang PetscFunctionReturn(0); 15191ad426b7SHong Zhang } 1520a4e85ca8SHong Zhang 1521e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 15223e97b6e8SHong Zhang { 15233e97b6e8SHong Zhang PetscErrorCode ierr; 15243e97b6e8SHong Zhang PetscInt j,M,N,row,col,ncols_u; 15253e97b6e8SHong Zhang const PetscInt *cols; 15263e97b6e8SHong Zhang PetscScalar zero=0.0; 15273e97b6e8SHong Zhang 15283e97b6e8SHong Zhang PetscFunctionBegin; 15293e97b6e8SHong Zhang ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 15303e97b6e8SHong 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); 15313e97b6e8SHong Zhang 15323e97b6e8SHong Zhang for (row=0; row<nrows; row++) { 15333e97b6e8SHong Zhang ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 15343e97b6e8SHong Zhang for (j=0; j<ncols_u; j++) { 15353e97b6e8SHong Zhang col = cols[j] + cstart; 15363e97b6e8SHong Zhang ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 15373e97b6e8SHong Zhang } 15383e97b6e8SHong Zhang ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 15393e97b6e8SHong Zhang } 15403e97b6e8SHong Zhang PetscFunctionReturn(0); 15413e97b6e8SHong Zhang } 15421ad426b7SHong Zhang 1543e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1544a4e85ca8SHong Zhang { 1545a4e85ca8SHong Zhang PetscErrorCode ierr; 1546f4431b8cSHong Zhang 1547a4e85ca8SHong Zhang PetscFunctionBegin; 1548a4e85ca8SHong Zhang if (Ju) { 1549a4e85ca8SHong Zhang ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1550a4e85ca8SHong Zhang } else { 1551a4e85ca8SHong Zhang ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1552a4e85ca8SHong Zhang } 1553a4e85ca8SHong Zhang PetscFunctionReturn(0); 1554a4e85ca8SHong Zhang } 1555a4e85ca8SHong Zhang 155624121865SAdrian 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. 155724121865SAdrian Maldonado */ 155824121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 155924121865SAdrian Maldonado { 156024121865SAdrian Maldonado PetscErrorCode ierr; 156124121865SAdrian Maldonado PetscInt i, size, dof; 156224121865SAdrian Maldonado PetscInt *glob2loc; 156324121865SAdrian Maldonado 156424121865SAdrian Maldonado PetscFunctionBegin; 156524121865SAdrian Maldonado ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr); 156624121865SAdrian Maldonado ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr); 156724121865SAdrian Maldonado 156824121865SAdrian Maldonado for (i = 0; i < size; i++) { 156924121865SAdrian Maldonado ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr); 157024121865SAdrian Maldonado dof = (dof >= 0) ? dof : -(dof + 1); 157124121865SAdrian Maldonado glob2loc[i] = dof; 157224121865SAdrian Maldonado } 157324121865SAdrian Maldonado 157424121865SAdrian Maldonado ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 157524121865SAdrian Maldonado #if 0 157624121865SAdrian Maldonado ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 157724121865SAdrian Maldonado #endif 157824121865SAdrian Maldonado PetscFunctionReturn(0); 157924121865SAdrian Maldonado } 158024121865SAdrian Maldonado 158101ad2aeeSHong Zhang #include <petsc/private/matimpl.h> 15821ad426b7SHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 15831ad426b7SHong Zhang { 15841ad426b7SHong Zhang PetscErrorCode ierr; 158524121865SAdrian Maldonado PetscMPIInt rank, size; 15861ad426b7SHong Zhang DM_Network *network = (DM_Network*) dm->data; 1587a4e85ca8SHong Zhang PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 1588840c2264SHong Zhang PetscInt cstart,ncols,j,e,v; 158924121865SAdrian Maldonado PetscBool ghost,ghost_vc,ghost2,isNest; 1590a4e85ca8SHong Zhang Mat Juser; 1591bfbc38dcSHong Zhang PetscSection sectionGlobal; 1592447d78afSSatish Balay PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 1593a4e85ca8SHong Zhang const PetscInt *edges,*cone; 15945cf7da58SHong Zhang MPI_Comm comm; 159524121865SAdrian Maldonado MatType mtype; 15965cf7da58SHong Zhang Vec vd_nz,vo_nz; 15975cf7da58SHong Zhang PetscInt *dnnz,*onnz; 15985cf7da58SHong Zhang PetscScalar *vdnz,*vonz; 15991ad426b7SHong Zhang 16001ad426b7SHong Zhang PetscFunctionBegin; 160124121865SAdrian Maldonado mtype = dm->mattype; 160224121865SAdrian Maldonado ierr = PetscStrcmp(mtype, MATNEST, &isNest);CHKERRQ(ierr); 160324121865SAdrian Maldonado 160424121865SAdrian Maldonado if (isNest) { 16050731d606SHong Zhang /* ierr = DMCreateMatrix_Network_Nest(); */ 160624121865SAdrian Maldonado PetscInt eDof, vDof; 160724121865SAdrian Maldonado Mat j11, j12, j21, j22, bA[2][2]; 160824121865SAdrian Maldonado ISLocalToGlobalMapping eISMap, vISMap; 160924121865SAdrian Maldonado 161024121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 161124121865SAdrian Maldonado ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 161224121865SAdrian Maldonado ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 161324121865SAdrian Maldonado 161424121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr); 161524121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr); 161624121865SAdrian Maldonado 161701ad2aeeSHong Zhang ierr = MatCreate(comm, &j11);CHKERRQ(ierr); 161824121865SAdrian Maldonado ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 161924121865SAdrian Maldonado ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr); 162024121865SAdrian Maldonado 162101ad2aeeSHong Zhang ierr = MatCreate(comm, &j12);CHKERRQ(ierr); 162224121865SAdrian Maldonado ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr); 162324121865SAdrian Maldonado ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr); 162424121865SAdrian Maldonado 162501ad2aeeSHong Zhang ierr = MatCreate(comm, &j21);CHKERRQ(ierr); 162624121865SAdrian Maldonado ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 162724121865SAdrian Maldonado ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr); 162824121865SAdrian Maldonado 162901ad2aeeSHong Zhang ierr = MatCreate(comm, &j22);CHKERRQ(ierr); 163024121865SAdrian Maldonado ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 163124121865SAdrian Maldonado ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr); 163224121865SAdrian Maldonado 16333f6a6bdaSHong Zhang bA[0][0] = j11; 16343f6a6bdaSHong Zhang bA[0][1] = j12; 16353f6a6bdaSHong Zhang bA[1][0] = j21; 16363f6a6bdaSHong Zhang bA[1][1] = j22; 163724121865SAdrian Maldonado 163824121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr); 163924121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr); 164024121865SAdrian Maldonado 164124121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr); 164224121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr); 164324121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr); 164424121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr); 164524121865SAdrian Maldonado 164624121865SAdrian Maldonado ierr = MatSetUp(j11);CHKERRQ(ierr); 164724121865SAdrian Maldonado ierr = MatSetUp(j12);CHKERRQ(ierr); 164824121865SAdrian Maldonado ierr = MatSetUp(j21);CHKERRQ(ierr); 164924121865SAdrian Maldonado ierr = MatSetUp(j22);CHKERRQ(ierr); 165024121865SAdrian Maldonado 165101ad2aeeSHong Zhang ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr); 165224121865SAdrian Maldonado ierr = MatSetUp(*J);CHKERRQ(ierr); 165324121865SAdrian Maldonado ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr); 165424121865SAdrian Maldonado ierr = MatDestroy(&j11);CHKERRQ(ierr); 165524121865SAdrian Maldonado ierr = MatDestroy(&j12);CHKERRQ(ierr); 165624121865SAdrian Maldonado ierr = MatDestroy(&j21);CHKERRQ(ierr); 165724121865SAdrian Maldonado ierr = MatDestroy(&j22);CHKERRQ(ierr); 165824121865SAdrian Maldonado 165924121865SAdrian Maldonado ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166024121865SAdrian Maldonado ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 166124121865SAdrian Maldonado 166224121865SAdrian Maldonado /* Free structures */ 166324121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr); 166424121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr); 166524121865SAdrian Maldonado 166624121865SAdrian Maldonado PetscFunctionReturn(0); 166724121865SAdrian Maldonado } else if (!network->userEdgeJacobian && !network->userVertexJacobian) { 1668a4e85ca8SHong Zhang /* user does not provide Jacobian blocks */ 1669bfbc38dcSHong Zhang ierr = DMCreateMatrix(network->plex,J);CHKERRQ(ierr); 1670bfbc38dcSHong Zhang ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 16711ad426b7SHong Zhang PetscFunctionReturn(0); 16721ad426b7SHong Zhang } 16731ad426b7SHong Zhang 1674bfbc38dcSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 16752a945128SHong Zhang ierr = DMGetDefaultGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 1676bfbc38dcSHong Zhang ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 1677bfbc38dcSHong Zhang ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 16782a945128SHong Zhang 16792a945128SHong Zhang ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 16802a945128SHong Zhang ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 168189898e50SHong Zhang 168289898e50SHong Zhang /* (1) Set matrix preallocation */ 168389898e50SHong Zhang /*------------------------------*/ 1684840c2264SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1685840c2264SHong Zhang ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr); 1686840c2264SHong Zhang ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr); 1687840c2264SHong Zhang ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr); 1688840c2264SHong Zhang ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr); 1689840c2264SHong Zhang ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr); 1690840c2264SHong Zhang 169189898e50SHong Zhang /* Set preallocation for edges */ 169289898e50SHong Zhang /*-----------------------------*/ 1693840c2264SHong Zhang ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 1694840c2264SHong Zhang 1695bdcb62a2SHong Zhang ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr); 1696840c2264SHong Zhang for (e=eStart; e<eEnd; e++) { 1697840c2264SHong Zhang /* Get row indices */ 1698840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 1699840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 1700840c2264SHong Zhang if (nrows) { 1701840c2264SHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 1702840c2264SHong Zhang 17035cf7da58SHong Zhang /* Set preallocation for conntected vertices */ 1704d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1705840c2264SHong Zhang for (v=0; v<2; v++) { 1706840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 1707840c2264SHong Zhang 17088675203cSHong Zhang if (network->Je) { 1709840c2264SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 17108675203cSHong Zhang } else Juser = NULL; 1711840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 17125cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1713840c2264SHong Zhang } 1714840c2264SHong Zhang 171589898e50SHong Zhang /* Set preallocation for edge self */ 1716840c2264SHong Zhang cstart = rstart; 17178675203cSHong Zhang if (network->Je) { 1718840c2264SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 17198675203cSHong Zhang } else Juser = NULL; 17205cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1721840c2264SHong Zhang } 1722840c2264SHong Zhang } 1723840c2264SHong Zhang 172489898e50SHong Zhang /* Set preallocation for vertices */ 172589898e50SHong Zhang /*--------------------------------*/ 1726840c2264SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 17278675203cSHong Zhang if (vEnd - vStart) vptr = network->Jvptr; 1728840c2264SHong Zhang 1729840c2264SHong Zhang for (v=vStart; v<vEnd; v++) { 1730840c2264SHong Zhang /* Get row indices */ 1731840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1732840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 1733840c2264SHong Zhang if (!nrows) continue; 1734840c2264SHong Zhang 1735bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1736bdcb62a2SHong Zhang if (ghost) { 1737bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1738bdcb62a2SHong Zhang } else { 1739bdcb62a2SHong Zhang rows_v = rows; 1740bdcb62a2SHong Zhang } 1741bdcb62a2SHong Zhang 1742bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1743840c2264SHong Zhang 1744840c2264SHong Zhang /* Get supporting edges and connected vertices */ 1745840c2264SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1746840c2264SHong Zhang 1747840c2264SHong Zhang for (e=0; e<nedges; e++) { 1748840c2264SHong Zhang /* Supporting edges */ 1749840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1750840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1751840c2264SHong Zhang 17528675203cSHong Zhang if (network->Jv) { 1753840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 17548675203cSHong Zhang } else Juser = NULL; 1755bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1756840c2264SHong Zhang 1757840c2264SHong Zhang /* Connected vertices */ 1758d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 1759840c2264SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 1760840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr); 1761840c2264SHong Zhang 1762840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1763840c2264SHong Zhang 17648675203cSHong Zhang if (network->Jv) { 1765840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 17668675203cSHong Zhang } else Juser = NULL; 1767e102a522SHong Zhang if (ghost_vc||ghost) { 1768e102a522SHong Zhang ghost2 = PETSC_TRUE; 1769e102a522SHong Zhang } else { 1770e102a522SHong Zhang ghost2 = PETSC_FALSE; 1771e102a522SHong Zhang } 1772e102a522SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr); 1773840c2264SHong Zhang } 1774840c2264SHong Zhang 177589898e50SHong Zhang /* Set preallocation for vertex self */ 1776840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1777840c2264SHong Zhang if (!ghost) { 1778840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 17798675203cSHong Zhang if (network->Jv) { 1780840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 17818675203cSHong Zhang } else Juser = NULL; 1782bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1783840c2264SHong Zhang } 1784bdcb62a2SHong Zhang if (ghost) { 1785bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1786bdcb62a2SHong Zhang } 1787840c2264SHong Zhang } 1788840c2264SHong Zhang 1789840c2264SHong Zhang ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr); 1790840c2264SHong Zhang ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr); 17915cf7da58SHong Zhang 17925cf7da58SHong Zhang ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr); 17935cf7da58SHong Zhang 17945cf7da58SHong Zhang ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr); 1795840c2264SHong Zhang ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr); 1796840c2264SHong Zhang 1797840c2264SHong Zhang ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr); 1798840c2264SHong Zhang ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr); 1799840c2264SHong Zhang for (j=0; j<localSize; j++) { 1800e102a522SHong Zhang dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 1801e102a522SHong Zhang onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 1802840c2264SHong Zhang } 1803840c2264SHong Zhang ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr); 1804840c2264SHong Zhang ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr); 1805840c2264SHong Zhang ierr = VecDestroy(&vd_nz);CHKERRQ(ierr); 1806840c2264SHong Zhang ierr = VecDestroy(&vo_nz);CHKERRQ(ierr); 1807840c2264SHong Zhang 18085cf7da58SHong Zhang ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr); 18095cf7da58SHong Zhang ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr); 18105cf7da58SHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 18115cf7da58SHong Zhang 18125cf7da58SHong Zhang ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr); 18135cf7da58SHong Zhang 181489898e50SHong Zhang /* (2) Set matrix entries for edges */ 181589898e50SHong Zhang /*----------------------------------*/ 18161ad426b7SHong Zhang for (e=eStart; e<eEnd; e++) { 1817bfbc38dcSHong Zhang /* Get row indices */ 18181ad426b7SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 181917df6e9eSHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 18204b976069SHong Zhang if (nrows) { 182117df6e9eSHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 18221ad426b7SHong Zhang 1823bfbc38dcSHong Zhang /* Set matrix entries for conntected vertices */ 1824d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1825bfbc38dcSHong Zhang for (v=0; v<2; v++) { 1826bfbc38dcSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 1827883e35e8SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 18283e97b6e8SHong Zhang 18298675203cSHong Zhang if (network->Je) { 1830a4e85ca8SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 18318675203cSHong Zhang } else Juser = NULL; 1832a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1833bfbc38dcSHong Zhang } 183417df6e9eSHong Zhang 1835bfbc38dcSHong Zhang /* Set matrix entries for edge self */ 18363e97b6e8SHong Zhang cstart = rstart; 18378675203cSHong Zhang if (network->Je) { 1838a4e85ca8SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 18398675203cSHong Zhang } else Juser = NULL; 1840a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 18411ad426b7SHong Zhang } 18424b976069SHong Zhang } 18431ad426b7SHong Zhang 1844bfbc38dcSHong Zhang /* Set matrix entries for vertices */ 184583b2e829SHong Zhang /*---------------------------------*/ 18461ad426b7SHong Zhang for (v=vStart; v<vEnd; v++) { 1847bfbc38dcSHong Zhang /* Get row indices */ 1848596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1849596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 18504b976069SHong Zhang if (!nrows) continue; 1851596e729fSHong Zhang 1852bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1853bdcb62a2SHong Zhang if (ghost) { 1854bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1855bdcb62a2SHong Zhang } else { 1856bdcb62a2SHong Zhang rows_v = rows; 1857bdcb62a2SHong Zhang } 1858bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1859596e729fSHong Zhang 1860bfbc38dcSHong Zhang /* Get supporting edges and connected vertices */ 1861596e729fSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1862596e729fSHong Zhang 1863596e729fSHong Zhang for (e=0; e<nedges; e++) { 1864bfbc38dcSHong Zhang /* Supporting edges */ 1865596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1866596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1867596e729fSHong Zhang 18688675203cSHong Zhang if (network->Jv) { 1869a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 18708675203cSHong Zhang } else Juser = NULL; 1871bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1872596e729fSHong Zhang 1873bfbc38dcSHong Zhang /* Connected vertices */ 1874d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 18752a945128SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 18762a945128SHong Zhang 187744aca652SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 187844aca652SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1879a4e85ca8SHong Zhang 18808675203cSHong Zhang if (network->Jv) { 1881a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 18828675203cSHong Zhang } else Juser = NULL; 1883bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1884596e729fSHong Zhang } 1885596e729fSHong Zhang 1886bfbc38dcSHong Zhang /* Set matrix entries for vertex self */ 18871ad426b7SHong Zhang if (!ghost) { 1888596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 18898675203cSHong Zhang if (network->Jv) { 1890a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 18918675203cSHong Zhang } else Juser = NULL; 1892bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr); 1893bdcb62a2SHong Zhang } 1894bdcb62a2SHong Zhang if (ghost) { 1895bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1896bdcb62a2SHong Zhang } 18971ad426b7SHong Zhang } 1898a4e85ca8SHong Zhang ierr = PetscFree(rows);CHKERRQ(ierr); 1899bdcb62a2SHong Zhang 19001ad426b7SHong Zhang ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19011ad426b7SHong Zhang ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1902dd6f46cdSHong Zhang 19035f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 19045f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19055f2c45f1SShri Abhyankar } 19065f2c45f1SShri Abhyankar 19075f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 19085f2c45f1SShri Abhyankar { 19095f2c45f1SShri Abhyankar PetscErrorCode ierr; 19105f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19112727e31bSShri Abhyankar PetscInt j; 19125f2c45f1SShri Abhyankar 19135f2c45f1SShri Abhyankar PetscFunctionBegin; 19148415c774SShri Abhyankar if (--network->refct > 0) PetscFunctionReturn(0); 191583b2e829SHong Zhang if (network->Je) { 191683b2e829SHong Zhang ierr = PetscFree(network->Je);CHKERRQ(ierr); 191783b2e829SHong Zhang } 191883b2e829SHong Zhang if (network->Jv) { 1919883e35e8SHong Zhang ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 192083b2e829SHong Zhang ierr = PetscFree(network->Jv);CHKERRQ(ierr); 19211ad426b7SHong Zhang } 192213c2a604SAdrian Maldonado 192313c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr); 192413c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr); 192513c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr); 192613c2a604SAdrian Maldonado if (network->vertex.sf) { 192713c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr); 192813c2a604SAdrian Maldonado } 192913c2a604SAdrian Maldonado /* edge */ 193013c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr); 193113c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr); 193213c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr); 193313c2a604SAdrian Maldonado if (network->edge.sf) { 193413c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr); 193513c2a604SAdrian Maldonado } 19365f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 19375f2c45f1SShri Abhyankar network->edges = NULL; 19385f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 19395f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 194083b2e829SHong Zhang 19412727e31bSShri Abhyankar for(j=0; j < network->nsubnet; j++) { 19422727e31bSShri Abhyankar ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr); 19432727e31bSShri Abhyankar ierr = PetscFree(network->subnet[j].vertices);CHKERRQ(ierr); 19442727e31bSShri Abhyankar } 1945e2aaf10cSShri Abhyankar ierr = PetscFree(network->subnet);CHKERRQ(ierr); 19465f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 19475f2c45f1SShri Abhyankar ierr = PetscFree(network->cvalue);CHKERRQ(ierr); 19485f2c45f1SShri Abhyankar ierr = PetscFree(network->header);CHKERRQ(ierr); 19495f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 19505f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19515f2c45f1SShri Abhyankar } 19525f2c45f1SShri Abhyankar 19535f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm, PetscViewer viewer) 19545f2c45f1SShri Abhyankar { 19555f2c45f1SShri Abhyankar PetscErrorCode ierr; 19565f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19575f2c45f1SShri Abhyankar 19585f2c45f1SShri Abhyankar PetscFunctionBegin; 19595f2c45f1SShri Abhyankar ierr = DMView(network->plex,viewer);CHKERRQ(ierr); 19605f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19615f2c45f1SShri Abhyankar } 19625f2c45f1SShri Abhyankar 19635f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 19645f2c45f1SShri Abhyankar { 19655f2c45f1SShri Abhyankar PetscErrorCode ierr; 19665f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19675f2c45f1SShri Abhyankar 19685f2c45f1SShri Abhyankar PetscFunctionBegin; 19695f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 19705f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19715f2c45f1SShri Abhyankar } 19725f2c45f1SShri Abhyankar 19735f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 19745f2c45f1SShri Abhyankar { 19755f2c45f1SShri Abhyankar PetscErrorCode ierr; 19765f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19775f2c45f1SShri Abhyankar 19785f2c45f1SShri Abhyankar PetscFunctionBegin; 19795f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 19805f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19815f2c45f1SShri Abhyankar } 19825f2c45f1SShri Abhyankar 19835f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 19845f2c45f1SShri Abhyankar { 19855f2c45f1SShri Abhyankar PetscErrorCode ierr; 19865f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19875f2c45f1SShri Abhyankar 19885f2c45f1SShri Abhyankar PetscFunctionBegin; 19895f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 19905f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19915f2c45f1SShri Abhyankar } 19925f2c45f1SShri Abhyankar 19935f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 19945f2c45f1SShri Abhyankar { 19955f2c45f1SShri Abhyankar PetscErrorCode ierr; 19965f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 19975f2c45f1SShri Abhyankar 19985f2c45f1SShri Abhyankar PetscFunctionBegin; 19995f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 20005f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20015f2c45f1SShri Abhyankar } 2002