1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 25f2c45f1SShri Abhyankar 35f2c45f1SShri Abhyankar /*@ 4556ed216SShri Abhyankar DMNetworkGetPlex - Gets the Plex DM associated with this network DM 5556ed216SShri Abhyankar 6556ed216SShri Abhyankar Not collective 7556ed216SShri Abhyankar 8556ed216SShri Abhyankar Input Parameters: 9556ed216SShri Abhyankar + netdm - the dm object 10556ed216SShri Abhyankar - plexmdm - the plex dm object 11556ed216SShri Abhyankar 12556ed216SShri Abhyankar Level: Advanced 13556ed216SShri Abhyankar 14556ed216SShri Abhyankar .seealso: DMNetworkCreate() 15556ed216SShri Abhyankar @*/ 16556ed216SShri Abhyankar PetscErrorCode DMNetworkGetPlex(DM netdm, DM *plexdm) 17556ed216SShri Abhyankar { 18556ed216SShri Abhyankar DM_Network *network = (DM_Network*) netdm->data; 19556ed216SShri Abhyankar 20556ed216SShri Abhyankar PetscFunctionBegin; 21556ed216SShri Abhyankar *plexdm = network->plex; 22556ed216SShri Abhyankar PetscFunctionReturn(0); 23556ed216SShri Abhyankar } 24556ed216SShri Abhyankar 25556ed216SShri Abhyankar /*@ 26e2aaf10cSShri Abhyankar DMNetworkSetSizes - Sets the number of subnetworks, local and global vertices and edges for each subnetwork. 275f2c45f1SShri Abhyankar 285f2c45f1SShri Abhyankar Collective on DM 295f2c45f1SShri Abhyankar 305f2c45f1SShri Abhyankar Input Parameters: 315f2c45f1SShri Abhyankar + dm - the dm object 32caf410d2SHong Zhang . Nsubnet - global number of subnetworks 33caf410d2SHong Zhang . nV - number of local vertices for each subnetwork 34caf410d2SHong Zhang . nE - number of local edges for each subnetwork 35caf410d2SHong Zhang . NsubnetCouple - global number of coupling subnetworks 36caf410d2SHong Zhang - nec - number of local edges for each coupling subnetwork 375f2c45f1SShri Abhyankar 38caf410d2SHong Zhang You cannot change the sizes once they have been set. nV, nE are arrays of length Nsubnet, and nec is array of length NsubnetCouple. 395f2c45f1SShri Abhyankar 401b266c99SBarry Smith Level: intermediate 411b266c99SBarry Smith 421b266c99SBarry Smith .seealso: DMNetworkCreate() 435f2c45f1SShri Abhyankar @*/ 44caf410d2SHong Zhang PetscErrorCode DMNetworkSetSizes(DM dm,PetscInt Nsubnet,PetscInt nV[], PetscInt nE[],PetscInt NsubnetCouple,PetscInt nec[]) 455f2c45f1SShri Abhyankar { 465f2c45f1SShri Abhyankar PetscErrorCode ierr; 475f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 48e2aaf10cSShri Abhyankar PetscInt a[2],b[2],i; 495f2c45f1SShri Abhyankar 505f2c45f1SShri Abhyankar PetscFunctionBegin; 515f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 52e2aaf10cSShri Abhyankar if (Nsubnet <= 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of subnetworks %D cannot be less than 1",Nsubnet); 537765340cSHong Zhang if (NsubnetCouple < 0) SETERRQ1(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of coupling subnetworks %D cannot be less than 0",NsubnetCouple); 547765340cSHong Zhang 55caf410d2SHong Zhang PetscValidLogicalCollectiveInt(dm,Nsubnet,2); 56caf410d2SHong Zhang if (NsubnetCouple) PetscValidLogicalCollectiveInt(dm,NsubnetCouple,5); 577765340cSHong Zhang if (network->nsubnet != 0) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network"); 587765340cSHong Zhang 59caf410d2SHong Zhang if (!nV || !nE) SETERRQ(PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"Local vertex size or edge size must be provided"); 60f025b11dSHong Zhang 61caf410d2SHong Zhang network->nsubnet = Nsubnet + NsubnetCouple; 62caf410d2SHong Zhang network->ncsubnet = NsubnetCouple; 63caf410d2SHong Zhang ierr = PetscCalloc1(Nsubnet+NsubnetCouple,&network->subnet);CHKERRQ(ierr); 64caf410d2SHong Zhang 65caf410d2SHong Zhang /* ---------------------------------------------------------- 66caf410d2SHong Zhang p=v or e; P=V or E 67caf410d2SHong Zhang subnet[0].pStart = 0 68caf410d2SHong Zhang subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i]) 69caf410d2SHong Zhang ----------------------------------------------------------------------- */ 70caf410d2SHong Zhang for (i=0; i < Nsubnet; i++) { 71caf410d2SHong Zhang /* Get global number of vertices and edges for subnet[i] */ 72caf410d2SHong Zhang a[0] = nV[i]; a[1] = nE[i]; /* local number of vertices (excluding ghost) and edges */ 737765340cSHong Zhang ierr = MPIU_Allreduce(a,b,2,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 747765340cSHong Zhang network->subnet[i].Nvtx = b[0]; network->subnet[i].Nedge = b[1]; 757765340cSHong Zhang 76caf410d2SHong Zhang network->subnet[i].nvtx = nV[i]; /* local nvtx, without ghost */ 777765340cSHong Zhang 78caf410d2SHong Zhang /* global subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */ 79caf410d2SHong Zhang network->subnet[i].vStart = network->NVertices; 807765340cSHong Zhang network->subnet[i].vEnd = network->subnet[i].vStart + network->subnet[i].Nvtx; 81caf410d2SHong Zhang 82caf410d2SHong Zhang network->nVertices += nV[i]; 837765340cSHong Zhang network->NVertices += network->subnet[i].Nvtx; 847765340cSHong Zhang 857765340cSHong Zhang network->subnet[i].nedge = nE[i]; 867765340cSHong Zhang network->subnet[i].eStart = network->nEdges; 87caf410d2SHong Zhang network->subnet[i].eEnd = network->subnet[i].eStart + nE[i]; 88caf410d2SHong Zhang network->nEdges += nE[i]; 89caf410d2SHong Zhang network->NEdges += network->subnet[i].Nedge; 90caf410d2SHong Zhang } 91caf410d2SHong Zhang 92caf410d2SHong Zhang /* coupling subnetwork */ 93caf410d2SHong Zhang for (; i < Nsubnet+NsubnetCouple; i++) { 94caf410d2SHong Zhang /* Get global number of coupling edges for subnet[i] */ 95caf410d2SHong Zhang ierr = MPIU_Allreduce(nec+(i-Nsubnet),&network->subnet[i].Nedge,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRQ(ierr); 96caf410d2SHong Zhang 97caf410d2SHong Zhang network->subnet[i].nvtx = 0; /* We design coupling subnetwork such that it does not have its own vertices */ 98caf410d2SHong Zhang network->subnet[i].vStart = network->nVertices; 99caf410d2SHong Zhang network->subnet[i].vEnd = network->subnet[i].vStart; 100caf410d2SHong Zhang 101caf410d2SHong Zhang network->subnet[i].nedge = nec[i-Nsubnet]; 102caf410d2SHong Zhang network->subnet[i].eStart = network->nEdges; 103caf410d2SHong Zhang network->subnet[i].eEnd = network->subnet[i].eStart + nec[i-Nsubnet]; 104caf410d2SHong Zhang network->nEdges += nec[i-Nsubnet]; 1057765340cSHong Zhang network->NEdges += network->subnet[i].Nedge; 1067765340cSHong Zhang } 1077765340cSHong Zhang PetscFunctionReturn(0); 1087765340cSHong Zhang } 1097765340cSHong Zhang 1105f2c45f1SShri Abhyankar /*@ 1115f2c45f1SShri Abhyankar DMNetworkSetEdgeList - Sets the list of local edges (vertex connectivity) for the network 1125f2c45f1SShri Abhyankar 1135f2c45f1SShri Abhyankar Logically collective on DM 1145f2c45f1SShri Abhyankar 1155f2c45f1SShri Abhyankar Input Parameters: 116e3e68989SHong Zhang + dm - the dm object 117e3e68989SHong Zhang . edgelist - list of edges for each subnetwork 118e3e68989SHong Zhang - edgelistCouple - list of edges for each coupling subnetwork 1195f2c45f1SShri Abhyankar 1205f2c45f1SShri Abhyankar Notes: 1215f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 1225f2c45f1SShri Abhyankar not be destroyed before the call to DMNetworkLayoutSetUp 1235f2c45f1SShri Abhyankar 1245f2c45f1SShri Abhyankar Level: intermediate 1255f2c45f1SShri Abhyankar 126e3e68989SHong Zhang Example usage: 127e3e68989SHong Zhang Consider the following 2 separate networks and a coupling network: 128e3e68989SHong Zhang 129e3e68989SHong Zhang .vb 130e3e68989SHong Zhang network 0: v0 -> v1 -> v2 -> v3 131e3e68989SHong Zhang network 1: v1 -> v2 -> v0 132e3e68989SHong Zhang coupling network: network 1: v2 -> network 0: v0 133e3e68989SHong Zhang .ve 134e3e68989SHong Zhang 135e3e68989SHong Zhang The resulting input 136e3e68989SHong Zhang edgelist[0] = [0 1 | 1 2 | 2 3]; 137e3e68989SHong Zhang edgelist[1] = [1 2 | 2 0] 138e3e68989SHong Zhang edgelistCouple[0] = [(network)1 (v)2 (network)0 (v)0]. 139e3e68989SHong Zhang 1405f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkSetSizes 1415f2c45f1SShri Abhyankar @*/ 1424e18019cSBarry Smith PetscErrorCode DMNetworkSetEdgeList(DM dm,PetscInt *edgelist[],PetscInt *edgelistCouple[]) 1435f2c45f1SShri Abhyankar { 1445f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 145e3e68989SHong Zhang PetscInt i,nsubnet,ncsubnet=network->ncsubnet; 1465f2c45f1SShri Abhyankar 1475f2c45f1SShri Abhyankar PetscFunctionBegin; 148e3e68989SHong Zhang nsubnet = network->nsubnet - ncsubnet; 149e3e68989SHong Zhang for(i=0; i < nsubnet; i++) { 150e2aaf10cSShri Abhyankar network->subnet[i].edgelist = edgelist[i]; 151e2aaf10cSShri Abhyankar } 152e3e68989SHong Zhang if (edgelistCouple) { 153e3e68989SHong Zhang PetscInt j; 154e3e68989SHong Zhang j = 0; 155e3e68989SHong Zhang nsubnet = network->nsubnet; 156e3e68989SHong Zhang while (i < nsubnet) network->subnet[i++].edgelist = edgelistCouple[j++]; 157e3e68989SHong Zhang } 1585f2c45f1SShri Abhyankar PetscFunctionReturn(0); 1595f2c45f1SShri Abhyankar } 1605f2c45f1SShri Abhyankar 1615f2c45f1SShri Abhyankar /*@ 1625f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 1635f2c45f1SShri Abhyankar 1645f2c45f1SShri Abhyankar Collective on DM 1655f2c45f1SShri Abhyankar 1665f2c45f1SShri Abhyankar Input Parameters 1675f2c45f1SShri Abhyankar . DM - the dmnetwork object 1685f2c45f1SShri Abhyankar 1695f2c45f1SShri Abhyankar Notes: 1705f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 1715f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 1725f2c45f1SShri Abhyankar 1735f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 1745f2c45f1SShri Abhyankar 1755f2c45f1SShri Abhyankar Level: intermediate 1765f2c45f1SShri Abhyankar 1775f2c45f1SShri Abhyankar .seealso: DMNetworkSetSizes, DMNetworkSetEdgeList 1785f2c45f1SShri Abhyankar @*/ 1795f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 1805f2c45f1SShri Abhyankar { 1815f2c45f1SShri Abhyankar PetscErrorCode ierr; 1825f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 183caf410d2SHong Zhang PetscInt numCorners=2,spacedim=2,dim = 1; /* One dimensional network */ 1843ebf9fb9SHong Zhang PetscReal *vertexcoords=NULL; 185caf410d2SHong Zhang PetscInt i,j,ctr,nsubnet,*eowners,np,*edges,*subnetvtx,vStart; 186caf410d2SHong Zhang PetscInt k,netid,vid, *vidxlTog,*edgelist_couple=NULL; 187caf410d2SHong Zhang const PetscInt *cone; 188caf410d2SHong Zhang MPI_Comm comm; 189caf410d2SHong Zhang PetscMPIInt size,rank; 1906500d4abSHong Zhang 1916500d4abSHong Zhang PetscFunctionBegin; 192caf410d2SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 193caf410d2SHong Zhang ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 194caf410d2SHong Zhang ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 1956500d4abSHong Zhang 196caf410d2SHong Zhang /* Create the local edgelist for the network by concatenating local input edgelists of the subnetworks */ 197caf410d2SHong Zhang ierr = PetscCalloc2(numCorners*network->nVertices,&vertexcoords,2*network->nEdges,&edges);CHKERRQ(ierr); 1987765340cSHong Zhang nsubnet = network->nsubnet - network->ncsubnet; 199caf410d2SHong Zhang ctr = 0; 2007765340cSHong Zhang for (i=0; i < nsubnet; i++) { 2016500d4abSHong Zhang for (j = 0; j < network->subnet[i].nedge; j++) { 202caf410d2SHong Zhang edges[2*ctr] = network->subnet[i].eStart + network->subnet[i].edgelist[2*j]; 203caf410d2SHong Zhang edges[2*ctr+1] = network->subnet[i].eStart + network->subnet[i].edgelist[2*j+1]; 2046500d4abSHong Zhang ctr++; 2056500d4abSHong Zhang } 2066500d4abSHong Zhang } 2077765340cSHong Zhang 208caf410d2SHong Zhang /* Append local coupling edgelists of the subnetworks */ 2097765340cSHong Zhang i = nsubnet; /* netid of coupling subnet */ 2107765340cSHong Zhang nsubnet = network->nsubnet; 2117765340cSHong Zhang while (i < nsubnet) { 212991cf414SHong Zhang edgelist_couple = network->subnet[i].edgelist; 2136500d4abSHong Zhang k = 0; 2146500d4abSHong Zhang for (j = 0; j < network->subnet[i].nedge; j++) { 2156500d4abSHong Zhang netid = edgelist_couple[k]; vid = edgelist_couple[k+1]; 216caf410d2SHong Zhang edges[2*ctr] = network->subnet[netid].vStart + vid; k += 2; 2176500d4abSHong Zhang 2186500d4abSHong Zhang netid = edgelist_couple[k]; vid = edgelist_couple[k+1]; 219caf410d2SHong Zhang edges[2*ctr+1] = network->subnet[netid].vStart + vid; k+=2; 2206500d4abSHong Zhang ctr++; 2216500d4abSHong Zhang } 2227765340cSHong Zhang i++; 2237765340cSHong Zhang } 224caf410d2SHong Zhang /* 225caf410d2SHong Zhang if (rank == 0) { 226caf410d2SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"[%d] edgelist:\n",rank); 2276500d4abSHong Zhang for(i=0; i < network->nEdges; i++) { 228caf410d2SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF,"[%D %D]",edges[2*i],edges[2*i+1]);CHKERRQ(ierr); 2296500d4abSHong Zhang printf("\n"); 2306500d4abSHong Zhang } 231caf410d2SHong Zhang } 232caf410d2SHong Zhang */ 2336500d4abSHong Zhang 234caf410d2SHong Zhang /* Create network->plex */ 235acdb140fSBarry Smith #if defined(PETSC_USE_64BIT_INDICES) 236acdb140fSBarry Smith { 237caf410d2SHong Zhang int *edges64; 238caf410d2SHong Zhang np = network->nEdges*numCorners; 239caf410d2SHong Zhang ierr = PetscMalloc1(np,&edges64);CHKERRQ(ierr); 240caf410d2SHong Zhang for (i=0; i<np; i++) edges64[i] = (int)edges[i]; 241caf410d2SHong Zhang 242caf410d2SHong Zhang if (size == 1) { 243caf410d2SHong Zhang ierr = DMPlexCreateFromCellList(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges64,spacedim,(const double*)vertexcoords,&network->plex);CHKERRQ(ierr); 244caf410d2SHong Zhang } else { 2453ebf9fb9SHong Zhang ierr = DMPlexCreateFromCellListParallel(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges64,spacedim,(const PetscReal*)vertexcoords,NULL,&network->plex);CHKERRQ(ierr); 246acdb140fSBarry Smith } 247caf410d2SHong Zhang ierr = PetscFree(edges64);CHKERRQ(ierr); 248acdb140fSBarry Smith } 249acdb140fSBarry Smith #else 250caf410d2SHong Zhang if (size == 1) { 251caf410d2SHong Zhang ierr = DMPlexCreateFromCellList(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges,spacedim,(const double*)vertexcoords,&network->plex);CHKERRQ(ierr); 252caf410d2SHong Zhang } else { 2533ebf9fb9SHong Zhang ierr = DMPlexCreateFromCellListParallel(comm,dim,network->nEdges,network->nVertices,numCorners,PETSC_FALSE,(const int*)edges,spacedim,(const PetscReal*)vertexcoords,NULL,&network->plex);CHKERRQ(ierr); 2546500d4abSHong Zhang } 255caf410d2SHong Zhang #endif 2566500d4abSHong Zhang 2576500d4abSHong Zhang ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 2586500d4abSHong Zhang ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 2596500d4abSHong Zhang ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 260caf410d2SHong Zhang vStart = network->vStart; 2616500d4abSHong Zhang 262caf410d2SHong Zhang ierr = PetscSectionCreate(comm,&network->DataSection);CHKERRQ(ierr); 263caf410d2SHong Zhang ierr = PetscSectionCreate(comm,&network->DofSection);CHKERRQ(ierr); 2646500d4abSHong Zhang ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 2656500d4abSHong Zhang ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 2666500d4abSHong Zhang 267caf410d2SHong Zhang network->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 268caf410d2SHong Zhang np = network->pEnd - network->pStart; 269caf410d2SHong Zhang ierr = PetscCalloc2(np,&network->header,np,&network->cvalue);CHKERRQ(ierr); 270caf410d2SHong Zhang 271caf410d2SHong Zhang /* Create vidxlTog: maps local vertex index to global index */ 272caf410d2SHong Zhang np = network->vEnd - vStart; 273caf410d2SHong Zhang ierr = PetscMalloc2(np,&vidxlTog,size+1,&eowners);CHKERRQ(ierr); 274caf410d2SHong Zhang ctr = 0; 275caf410d2SHong Zhang for (i=network->eStart; i<network->eEnd; i++) { 276caf410d2SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,i,&cone);CHKERRQ(ierr); 277caf410d2SHong Zhang vidxlTog[cone[0] - vStart] = edges[2*ctr]; 278caf410d2SHong Zhang vidxlTog[cone[1] - vStart] = edges[2*ctr+1]; 279caf410d2SHong Zhang ctr++; 280caf410d2SHong Zhang } 281caf410d2SHong Zhang ierr = PetscFree2(vertexcoords,edges);CHKERRQ(ierr); 282caf410d2SHong Zhang 2836500d4abSHong Zhang /* Create vertices and edges array for the subnetworks */ 2846500d4abSHong Zhang for (j=0; j < network->nsubnet; j++) { 2856500d4abSHong Zhang ierr = PetscCalloc1(network->subnet[j].nedge,&network->subnet[j].edges);CHKERRQ(ierr); 286caf410d2SHong Zhang 2876500d4abSHong Zhang /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 2886500d4abSHong Zhang These get updated when the vertices and edges are added. */ 289caf410d2SHong Zhang network->subnet[j].nvtx = 0; 290caf410d2SHong Zhang network->subnet[j].nedge = 0; 2916500d4abSHong Zhang } 292caf410d2SHong Zhang ierr = PetscCalloc1(np,&network->subnetvtx);CHKERRQ(ierr); 2936500d4abSHong Zhang 294caf410d2SHong Zhang 295caf410d2SHong Zhang /* Get edge ownership */ 296caf410d2SHong Zhang np = network->eEnd - network->eStart; 297caf410d2SHong Zhang ierr = MPI_Allgather(&np,1,MPIU_INT,eowners+1,1,MPIU_INT,comm);CHKERRQ(ierr); 298caf410d2SHong Zhang eowners[0] = 0; 299caf410d2SHong Zhang for (i=2; i<=size; i++) eowners[i] += eowners[i-1]; 300caf410d2SHong Zhang 301caf410d2SHong Zhang i = 0; j = 0; 302caf410d2SHong Zhang while (i < np) { /* local edges, including coupling edges */ 303caf410d2SHong Zhang network->header[i].index = i + eowners[rank]; /* Global edge index */ 304caf410d2SHong Zhang 305caf410d2SHong Zhang if (j < network->nsubnet && i < network->subnet[j].eEnd) { 3066500d4abSHong Zhang network->header[i].subnetid = j; /* Subnetwork id */ 3076500d4abSHong Zhang network->subnet[j].edges[network->subnet[j].nedge++] = i; 308caf410d2SHong Zhang 309caf410d2SHong Zhang network->header[i].ndata = 0; 310caf410d2SHong Zhang ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 311caf410d2SHong Zhang network->header[i].offset[0] = 0; 312caf410d2SHong Zhang i++; 313caf410d2SHong Zhang } 314caf410d2SHong Zhang if (i >= network->subnet[j].eEnd) j++; 315caf410d2SHong Zhang } 316caf410d2SHong Zhang 317caf410d2SHong Zhang /* Count network->subnet[*].nvtx */ 318caf410d2SHong Zhang for (i=vStart; i<network->vEnd; i++) { /* local vertices, including ghosts */ 319caf410d2SHong Zhang k = vidxlTog[i-vStart]; 320caf410d2SHong Zhang for (j=0; j < network->nsubnet; j++) { 321caf410d2SHong Zhang if (network->subnet[j].vStart <= k && k < network->subnet[j].vEnd) { 322caf410d2SHong Zhang network->subnet[j].nvtx++; 3236500d4abSHong Zhang break; 3246500d4abSHong Zhang } 3256500d4abSHong Zhang } 3266500d4abSHong Zhang } 3276500d4abSHong Zhang 328caf410d2SHong Zhang /* Set network->subnet[*].vertices on array network->subnetvtx */ 329caf410d2SHong Zhang subnetvtx = network->subnetvtx; 3306500d4abSHong Zhang for (j=0; j<network->nsubnet; j++) { 331caf410d2SHong Zhang network->subnet[j].vertices = subnetvtx; 332caf410d2SHong Zhang subnetvtx += network->subnet[j].nvtx; 333caf410d2SHong Zhang network->subnet[j].nvtx = 0; 334caf410d2SHong Zhang } 335caf410d2SHong Zhang 336caf410d2SHong Zhang /* Set vertex array for the subnetworks */ 337caf410d2SHong Zhang for (i=vStart; i<network->vEnd; i++) { /* local vertices, including ghosts */ 338caf410d2SHong Zhang network->header[i].index = vidxlTog[i-vStart]; /* Global vertex index */ 339caf410d2SHong Zhang 340caf410d2SHong Zhang k = vidxlTog[i-vStart]; 341caf410d2SHong Zhang for (j=0; j < network->nsubnet; j++) { 342caf410d2SHong Zhang if (network->subnet[j].vStart <= k && k < network->subnet[j].vEnd) { 3436500d4abSHong Zhang network->header[i].subnetid = j; 3446500d4abSHong Zhang network->subnet[j].vertices[network->subnet[j].nvtx++] = i; 3456500d4abSHong Zhang break; 3466500d4abSHong Zhang } 3476500d4abSHong Zhang } 348caf410d2SHong Zhang 3496500d4abSHong Zhang network->header[i].ndata = 0; 3506500d4abSHong Zhang ierr = PetscSectionAddDof(network->DataSection,i,network->dataheadersize);CHKERRQ(ierr); 351caf410d2SHong Zhang network->header[i].offset[0] = 0; 3526500d4abSHong Zhang } 3536500d4abSHong Zhang 354caf410d2SHong Zhang ierr = PetscFree2(vidxlTog,eowners);CHKERRQ(ierr); 3555f2c45f1SShri Abhyankar PetscFunctionReturn(0); 3565f2c45f1SShri Abhyankar } 3575f2c45f1SShri Abhyankar 35894ef8ddeSSatish Balay /*@C 3592727e31bSShri Abhyankar DMNetworkGetSubnetworkInfo - Returns the info for the subnetwork 3602727e31bSShri Abhyankar 3612727e31bSShri Abhyankar Input Parameters 362caf410d2SHong Zhang + dm - the DM object 3632727e31bSShri Abhyankar - id - the ID (integer) of the subnetwork 3642727e31bSShri Abhyankar 3652727e31bSShri Abhyankar Output Parameters 3662727e31bSShri Abhyankar + nv - number of vertices (local) 3672727e31bSShri Abhyankar . ne - number of edges (local) 3682727e31bSShri Abhyankar . vtx - local vertices for this subnetwork 3692727e31bSShri Abhyankar . edge - local edges for this subnetwork 3702727e31bSShri Abhyankar 3712727e31bSShri Abhyankar Notes: 3722727e31bSShri Abhyankar Cannot call this routine before DMNetworkLayoutSetup() 3732727e31bSShri Abhyankar 37406dd6b0eSSatish Balay Level: intermediate 37506dd6b0eSSatish Balay 3762727e31bSShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 3772727e31bSShri Abhyankar @*/ 378caf410d2SHong Zhang PetscErrorCode DMNetworkGetSubnetworkInfo(DM dm,PetscInt id,PetscInt *nv, PetscInt *ne,const PetscInt **vtx, const PetscInt **edge) 3792727e31bSShri Abhyankar { 380caf410d2SHong Zhang DM_Network *network = (DM_Network*)dm->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 391caf410d2SHong Zhang DMNetworkGetSubnetworkCoupleInfo - Returns the info for the coupling subnetwork 392caf410d2SHong Zhang 393caf410d2SHong Zhang Input Parameters 394caf410d2SHong Zhang + dm - the DM object 395caf410d2SHong Zhang - id - the ID (integer) of the coupling subnetwork 396caf410d2SHong Zhang 397caf410d2SHong Zhang Output Parameters 398caf410d2SHong Zhang + ne - number of edges (local) 399caf410d2SHong Zhang - edge - local edges for this coupling subnetwork 400caf410d2SHong Zhang 401caf410d2SHong Zhang Notes: 402caf410d2SHong Zhang Cannot call this routine before DMNetworkLayoutSetup() 403caf410d2SHong Zhang 404caf410d2SHong Zhang Level: intermediate 405caf410d2SHong Zhang 406caf410d2SHong Zhang .seealso: DMNetworkGetSubnetworkInfo, DMNetworkLayoutSetUp, DMNetworkCreate 407caf410d2SHong Zhang @*/ 408caf410d2SHong Zhang PetscErrorCode DMNetworkGetSubnetworkCoupleInfo(DM dm,PetscInt id,PetscInt *ne,const PetscInt **edge) 409caf410d2SHong Zhang { 410caf410d2SHong Zhang DM_Network *net = (DM_Network*)dm->data; 411caf410d2SHong Zhang PetscInt id1 = id + net->nsubnet - net->ncsubnet; 412caf410d2SHong Zhang 413caf410d2SHong Zhang PetscFunctionBegin; 414caf410d2SHong Zhang *ne = net->subnet[id1].nedge; 415caf410d2SHong Zhang *edge = net->subnet[id1].edges; 416caf410d2SHong Zhang PetscFunctionReturn(0); 417caf410d2SHong Zhang } 418caf410d2SHong Zhang 419caf410d2SHong Zhang /*@C 4205f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 4215f2c45f1SShri Abhyankar 4225f2c45f1SShri Abhyankar Logically collective on DM 4235f2c45f1SShri Abhyankar 4245f2c45f1SShri Abhyankar Input Parameters 4255f2c45f1SShri Abhyankar + dm - the network object 4265f2c45f1SShri Abhyankar . name - the component name 4275f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 4285f2c45f1SShri Abhyankar 4295f2c45f1SShri Abhyankar Output Parameters 4305f2c45f1SShri Abhyankar . key - an integer key that defines the component 4315f2c45f1SShri Abhyankar 4325f2c45f1SShri Abhyankar Notes 4335f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 4345f2c45f1SShri Abhyankar 4355f2c45f1SShri Abhyankar Level: intermediate 4365f2c45f1SShri Abhyankar 4375f2c45f1SShri Abhyankar .seealso: DMNetworkLayoutSetUp, DMNetworkCreate 4385f2c45f1SShri Abhyankar @*/ 439caf410d2SHong Zhang PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,size_t size,PetscInt *key) 4405f2c45f1SShri Abhyankar { 4415f2c45f1SShri Abhyankar PetscErrorCode ierr; 4425f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 4435f2c45f1SShri Abhyankar DMNetworkComponent *component=&network->component[network->ncomponent]; 4445f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 4455f2c45f1SShri Abhyankar PetscInt i; 4465f2c45f1SShri Abhyankar 4475f2c45f1SShri Abhyankar PetscFunctionBegin; 4485f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 4495f2c45f1SShri Abhyankar ierr = PetscStrcmp(component->name,name,&flg);CHKERRQ(ierr); 4505f2c45f1SShri Abhyankar if (flg) { 4515f2c45f1SShri Abhyankar *key = i; 4525f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4535f2c45f1SShri Abhyankar } 4546d64e262SShri Abhyankar } 4556d64e262SShri Abhyankar if(network->ncomponent == MAX_COMPONENTS) { 4566d64e262SShri Abhyankar SETERRQ1(PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_OUTOFRANGE,"Number of components registered exceeds the max %D",MAX_COMPONENTS); 4575f2c45f1SShri Abhyankar } 4585f2c45f1SShri Abhyankar 4595f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 4605f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 4615f2c45f1SShri Abhyankar *key = network->ncomponent; 4625f2c45f1SShri Abhyankar network->ncomponent++; 4635f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4645f2c45f1SShri Abhyankar } 4655f2c45f1SShri Abhyankar 4665f2c45f1SShri Abhyankar /*@ 4675f2c45f1SShri Abhyankar DMNetworkGetVertexRange - Get the bounds [start, end) for the vertices. 4685f2c45f1SShri Abhyankar 4695f2c45f1SShri Abhyankar Not Collective 4705f2c45f1SShri Abhyankar 4715f2c45f1SShri Abhyankar Input Parameters: 4725f2c45f1SShri Abhyankar + dm - The DMNetwork object 4735f2c45f1SShri Abhyankar 4745f2c45f1SShri Abhyankar Output Paramters: 4755f2c45f1SShri Abhyankar + vStart - The first vertex point 4765f2c45f1SShri Abhyankar - vEnd - One beyond the last vertex point 4775f2c45f1SShri Abhyankar 4785f2c45f1SShri Abhyankar Level: intermediate 4795f2c45f1SShri Abhyankar 4805f2c45f1SShri Abhyankar .seealso: DMNetworkGetEdgeRange 4815f2c45f1SShri Abhyankar @*/ 4825f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 4835f2c45f1SShri Abhyankar { 4845f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 4855f2c45f1SShri Abhyankar 4865f2c45f1SShri Abhyankar PetscFunctionBegin; 4875f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 4885f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 4895f2c45f1SShri Abhyankar PetscFunctionReturn(0); 4905f2c45f1SShri Abhyankar } 4915f2c45f1SShri Abhyankar 4925f2c45f1SShri Abhyankar /*@ 4935f2c45f1SShri Abhyankar DMNetworkGetEdgeRange - Get the bounds [start, end) for the edges. 4945f2c45f1SShri Abhyankar 4955f2c45f1SShri Abhyankar Not Collective 4965f2c45f1SShri Abhyankar 4975f2c45f1SShri Abhyankar Input Parameters: 4985f2c45f1SShri Abhyankar + dm - The DMNetwork object 4995f2c45f1SShri Abhyankar 5005f2c45f1SShri Abhyankar Output Paramters: 5015f2c45f1SShri Abhyankar + eStart - The first edge point 5025f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 5035f2c45f1SShri Abhyankar 5045f2c45f1SShri Abhyankar Level: intermediate 5055f2c45f1SShri Abhyankar 5065f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange 5075f2c45f1SShri Abhyankar @*/ 5085f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 5095f2c45f1SShri Abhyankar { 5105f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 5115f2c45f1SShri Abhyankar 5125f2c45f1SShri Abhyankar PetscFunctionBegin; 5135f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 5145f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 5155f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5165f2c45f1SShri Abhyankar } 5175f2c45f1SShri Abhyankar 5187b6afd5bSHong Zhang /*@ 519e85e6aecSHong Zhang DMNetworkGetGlobalEdgeIndex - Get the user global numbering for the edge. 5207b6afd5bSHong Zhang 5217b6afd5bSHong Zhang Not Collective 5227b6afd5bSHong Zhang 5237b6afd5bSHong Zhang Input Parameters: 5247b6afd5bSHong Zhang + dm - DMNetwork object 525e85e6aecSHong Zhang - p - edge point 5267b6afd5bSHong Zhang 5277b6afd5bSHong Zhang Output Paramters: 528e85e6aecSHong Zhang . index - user global numbering for the edge 5297b6afd5bSHong Zhang 5307b6afd5bSHong Zhang Level: intermediate 5317b6afd5bSHong Zhang 532e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalVertexIndex 5337b6afd5bSHong Zhang @*/ 534e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index) 5357b6afd5bSHong Zhang { 5367b6afd5bSHong Zhang PetscErrorCode ierr; 5377b6afd5bSHong Zhang DM_Network *network = (DM_Network*)dm->data; 5387b6afd5bSHong Zhang PetscInt offsetp; 5397b6afd5bSHong Zhang DMNetworkComponentHeader header; 5407b6afd5bSHong Zhang 5417b6afd5bSHong Zhang PetscFunctionBegin; 542caf410d2SHong Zhang if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first"); 5437b6afd5bSHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 5447b6afd5bSHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 545e85e6aecSHong Zhang *index = header->index; 5467b6afd5bSHong Zhang PetscFunctionReturn(0); 5477b6afd5bSHong Zhang } 5487b6afd5bSHong Zhang 5495f2c45f1SShri Abhyankar /*@ 550e85e6aecSHong Zhang DMNetworkGetGlobalVertexIndex - Get the user global numbering for the vertex. 551e85e6aecSHong Zhang 552e85e6aecSHong Zhang Not Collective 553e85e6aecSHong Zhang 554e85e6aecSHong Zhang Input Parameters: 555e85e6aecSHong Zhang + dm - DMNetwork object 556e85e6aecSHong Zhang - p - vertex point 557e85e6aecSHong Zhang 558e85e6aecSHong Zhang Output Paramters: 559e85e6aecSHong Zhang . index - user global numbering for the vertex 560e85e6aecSHong Zhang 561e85e6aecSHong Zhang Level: intermediate 562e85e6aecSHong Zhang 563e85e6aecSHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex 564e85e6aecSHong Zhang @*/ 565e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index) 566e85e6aecSHong Zhang { 567e85e6aecSHong Zhang PetscErrorCode ierr; 568e85e6aecSHong Zhang DM_Network *network = (DM_Network*)dm->data; 569e85e6aecSHong Zhang PetscInt offsetp; 570e85e6aecSHong Zhang DMNetworkComponentHeader header; 571e85e6aecSHong Zhang 572e85e6aecSHong Zhang PetscFunctionBegin; 573caf410d2SHong Zhang if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first"); 574e85e6aecSHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 575e85e6aecSHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 576e85e6aecSHong Zhang *index = header->index; 577e85e6aecSHong Zhang PetscFunctionReturn(0); 578e85e6aecSHong Zhang } 579e85e6aecSHong Zhang 580c3b11c7cSShri Abhyankar /* 581c3b11c7cSShri Abhyankar DMNetworkGetComponentKeyOffset - Gets the type along with the offset for indexing the 582c3b11c7cSShri Abhyankar component value from the component data array 583c3b11c7cSShri Abhyankar 584c3b11c7cSShri Abhyankar Not Collective 585c3b11c7cSShri Abhyankar 586c3b11c7cSShri Abhyankar Input Parameters: 587c3b11c7cSShri Abhyankar + dm - The DMNetwork object 588c3b11c7cSShri Abhyankar . p - vertex/edge point 589c3b11c7cSShri Abhyankar - compnum - component number 590c3b11c7cSShri Abhyankar 591c3b11c7cSShri Abhyankar Output Parameters: 592c3b11c7cSShri Abhyankar + compkey - the key obtained when registering the component 593c3b11c7cSShri Abhyankar - offset - offset into the component data array associated with the vertex/edge point 594c3b11c7cSShri Abhyankar 595c3b11c7cSShri Abhyankar Notes: 596c3b11c7cSShri Abhyankar Typical usage: 597c3b11c7cSShri Abhyankar 598c3b11c7cSShri Abhyankar DMNetworkGetComponentDataArray(dm, &arr); 599c3b11c7cSShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 600c3b11c7cSShri Abhyankar Loop over vertices or edges 601c3b11c7cSShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 602c3b11c7cSShri Abhyankar Loop over numcomps 603c3b11c7cSShri Abhyankar DMNetworkGetComponentKeyOffset(dm,v,compnum,&key,&offset); 604c3b11c7cSShri Abhyankar compdata = (UserCompDataType)(arr+offset); 605c3b11c7cSShri Abhyankar 606c3b11c7cSShri Abhyankar Level: intermediate 607c3b11c7cSShri Abhyankar 608c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetComponentDataArray, 609c3b11c7cSShri Abhyankar */ 610c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponentKeyOffset(DM dm,PetscInt p, PetscInt compnum, PetscInt *compkey, PetscInt *offset) 611c3b11c7cSShri Abhyankar { 612c3b11c7cSShri Abhyankar PetscErrorCode ierr; 613c3b11c7cSShri Abhyankar PetscInt offsetp; 614c3b11c7cSShri Abhyankar DMNetworkComponentHeader header; 615c3b11c7cSShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 616c3b11c7cSShri Abhyankar 617c3b11c7cSShri Abhyankar PetscFunctionBegin; 618c3b11c7cSShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 619c3b11c7cSShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 620c3b11c7cSShri Abhyankar if (compkey) *compkey = header->key[compnum]; 621c3b11c7cSShri Abhyankar if (offset) *offset = offsetp+network->dataheadersize+header->offset[compnum]; 622c3b11c7cSShri Abhyankar PetscFunctionReturn(0); 623c3b11c7cSShri Abhyankar } 624c3b11c7cSShri Abhyankar 625c3b11c7cSShri Abhyankar /*@ 626c3b11c7cSShri Abhyankar DMNetworkGetComponent - Returns the network component and its key 627c3b11c7cSShri Abhyankar 628c3b11c7cSShri Abhyankar Not Collective 629c3b11c7cSShri Abhyankar 630c3b11c7cSShri Abhyankar Input Parameters 631c3b11c7cSShri Abhyankar + dm - DMNetwork object 632c3b11c7cSShri Abhyankar . p - edge or vertex point 633c3b11c7cSShri Abhyankar - compnum - component number 634c3b11c7cSShri Abhyankar 635c3b11c7cSShri Abhyankar Output Parameters: 636c3b11c7cSShri Abhyankar + compkey - the key set for this computing during registration 637c3b11c7cSShri Abhyankar - component - the component data 638c3b11c7cSShri Abhyankar 639c3b11c7cSShri Abhyankar Notes: 640c3b11c7cSShri Abhyankar Typical usage: 641c3b11c7cSShri Abhyankar 642c3b11c7cSShri Abhyankar DMNetworkGetVertex/EdgeRange(dm,&Start,&End); 643c3b11c7cSShri Abhyankar Loop over vertices or edges 644c3b11c7cSShri Abhyankar DMNetworkGetNumComponents(dm,v,&numcomps); 645c3b11c7cSShri Abhyankar Loop over numcomps 646c3b11c7cSShri Abhyankar DMNetworkGetComponent(dm,v,compnum,&key,&component); 647c3b11c7cSShri Abhyankar 648c3b11c7cSShri Abhyankar Level: intermediate 649c3b11c7cSShri Abhyankar 650c3b11c7cSShri Abhyankar .seealso: DMNetworkGetNumComponents, DMNetworkGetVariableOffset 651c3b11c7cSShri Abhyankar @*/ 652c3b11c7cSShri Abhyankar PetscErrorCode DMNetworkGetComponent(DM dm, PetscInt p, PetscInt compnum, PetscInt *key, void **component) 653c3b11c7cSShri Abhyankar { 654c3b11c7cSShri Abhyankar PetscErrorCode ierr; 655c3b11c7cSShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 656e108cb99SStefano Zampini PetscInt offsetd = 0; 657c3b11c7cSShri Abhyankar 658c3b11c7cSShri Abhyankar PetscFunctionBegin; 659c3b11c7cSShri Abhyankar ierr = DMNetworkGetComponentKeyOffset(dm,p,compnum,key,&offsetd);CHKERRQ(ierr); 660c3b11c7cSShri Abhyankar *component = network->componentdataarray+offsetd; 661c3b11c7cSShri Abhyankar PetscFunctionReturn(0); 662c3b11c7cSShri Abhyankar } 663c3b11c7cSShri Abhyankar 664e85e6aecSHong Zhang /*@ 665325661f6SSatish Balay DMNetworkAddComponent - Adds a network component at the given point (vertex/edge) 6665f2c45f1SShri Abhyankar 6675f2c45f1SShri Abhyankar Not Collective 6685f2c45f1SShri Abhyankar 6695f2c45f1SShri Abhyankar Input Parameters: 6705f2c45f1SShri Abhyankar + dm - The DMNetwork object 6715f2c45f1SShri Abhyankar . p - vertex/edge point 6725f2c45f1SShri Abhyankar . componentkey - component key returned while registering the component 6735f2c45f1SShri Abhyankar - compvalue - pointer to the data structure for the component 6745f2c45f1SShri Abhyankar 6755f2c45f1SShri Abhyankar Level: intermediate 6765f2c45f1SShri Abhyankar 6775f2c45f1SShri Abhyankar .seealso: DMNetworkGetVertexRange, DMNetworkGetEdgeRange, DMNetworkRegisterComponent 6785f2c45f1SShri Abhyankar @*/ 6795f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddComponent(DM dm, PetscInt p,PetscInt componentkey,void* compvalue) 6805f2c45f1SShri Abhyankar { 6815f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 68243a39a44SBarry Smith DMNetworkComponent *component = &network->component[componentkey]; 6835f2c45f1SShri Abhyankar DMNetworkComponentHeader header = &network->header[p]; 6845f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue = &network->cvalue[p]; 6855f2c45f1SShri Abhyankar PetscErrorCode ierr; 6865f2c45f1SShri Abhyankar 6875f2c45f1SShri Abhyankar PetscFunctionBegin; 688fa58f0a9SHong 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); 689fa58f0a9SHong Zhang 69043a39a44SBarry Smith header->size[header->ndata] = component->size; 69143a39a44SBarry Smith ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 6925f2c45f1SShri Abhyankar header->key[header->ndata] = componentkey; 6935f2c45f1SShri Abhyankar if (header->ndata != 0) header->offset[header->ndata] = header->offset[header->ndata-1] + header->size[header->ndata-1]; 6945f2c45f1SShri Abhyankar 6955f2c45f1SShri Abhyankar cvalue->data[header->ndata] = (void*)compvalue; 6965f2c45f1SShri Abhyankar header->ndata++; 6975f2c45f1SShri Abhyankar PetscFunctionReturn(0); 6985f2c45f1SShri Abhyankar } 6995f2c45f1SShri Abhyankar 7005f2c45f1SShri Abhyankar /*@ 7015f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 7025f2c45f1SShri Abhyankar 7035f2c45f1SShri Abhyankar Not Collective 7045f2c45f1SShri Abhyankar 7055f2c45f1SShri Abhyankar Input Parameters: 7065f2c45f1SShri Abhyankar + dm - The DMNetwork object 7075f2c45f1SShri Abhyankar . p - vertex/edge point 7085f2c45f1SShri Abhyankar 7095f2c45f1SShri Abhyankar Output Parameters: 7105f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 7115f2c45f1SShri Abhyankar 7125f2c45f1SShri Abhyankar Level: intermediate 7135f2c45f1SShri Abhyankar 7145f2c45f1SShri Abhyankar .seealso: DMNetworkRegisterComponent, DMNetworkAddComponent 7155f2c45f1SShri Abhyankar @*/ 7165f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 7175f2c45f1SShri Abhyankar { 7185f2c45f1SShri Abhyankar PetscErrorCode ierr; 7195f2c45f1SShri Abhyankar PetscInt offset; 7205f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7215f2c45f1SShri Abhyankar 7225f2c45f1SShri Abhyankar PetscFunctionBegin; 7235f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 7245f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 7255f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7265f2c45f1SShri Abhyankar } 7275f2c45f1SShri Abhyankar 7285f2c45f1SShri Abhyankar /*@ 7295f2c45f1SShri Abhyankar DMNetworkGetVariableOffset - Get the offset for accessing the variable associated with the given vertex/edge from the local vector. 7305f2c45f1SShri Abhyankar 7315f2c45f1SShri Abhyankar Not Collective 7325f2c45f1SShri Abhyankar 7335f2c45f1SShri Abhyankar Input Parameters: 7345f2c45f1SShri Abhyankar + dm - The DMNetwork object 7355f2c45f1SShri Abhyankar - p - the edge/vertex point 7365f2c45f1SShri Abhyankar 7375f2c45f1SShri Abhyankar Output Parameters: 7385f2c45f1SShri Abhyankar . offset - the offset 7395f2c45f1SShri Abhyankar 7405f2c45f1SShri Abhyankar Level: intermediate 7415f2c45f1SShri Abhyankar 7425f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 7435f2c45f1SShri Abhyankar @*/ 7445f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableOffset(DM dm,PetscInt p,PetscInt *offset) 7455f2c45f1SShri Abhyankar { 7465f2c45f1SShri Abhyankar PetscErrorCode ierr; 7475f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7485f2c45f1SShri Abhyankar 7495f2c45f1SShri Abhyankar PetscFunctionBegin; 7505f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultSection,p,offset);CHKERRQ(ierr); 7515f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7525f2c45f1SShri Abhyankar } 7535f2c45f1SShri Abhyankar 7545f2c45f1SShri Abhyankar /*@ 7555f2c45f1SShri Abhyankar DMNetworkGetVariableGlobalOffset - Get the global offset for the variable associated with the given vertex/edge from the global vector. 7565f2c45f1SShri Abhyankar 7575f2c45f1SShri Abhyankar Not Collective 7585f2c45f1SShri Abhyankar 7595f2c45f1SShri Abhyankar Input Parameters: 7605f2c45f1SShri Abhyankar + dm - The DMNetwork object 7615f2c45f1SShri Abhyankar - p - the edge/vertex point 7625f2c45f1SShri Abhyankar 7635f2c45f1SShri Abhyankar Output Parameters: 7645f2c45f1SShri Abhyankar . offsetg - the offset 7655f2c45f1SShri Abhyankar 7665f2c45f1SShri Abhyankar Level: intermediate 7675f2c45f1SShri Abhyankar 7685f2c45f1SShri Abhyankar .seealso: DMNetworkGetVariableOffset, DMGetLocalVector 7695f2c45f1SShri Abhyankar @*/ 7705f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVariableGlobalOffset(DM dm,PetscInt p,PetscInt *offsetg) 7715f2c45f1SShri Abhyankar { 7725f2c45f1SShri Abhyankar PetscErrorCode ierr; 7735f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 7745f2c45f1SShri Abhyankar 7755f2c45f1SShri Abhyankar PetscFunctionBegin; 7765f78ed8bSShri Abhyankar ierr = PetscSectionGetOffset(network->plex->defaultGlobalSection,p,offsetg);CHKERRQ(ierr); 7776fefedf4SHong Zhang if (*offsetg < 0) *offsetg = -(*offsetg + 1); /* Convert to actual global offset for ghost vertex */ 7785f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7795f2c45f1SShri Abhyankar } 7805f2c45f1SShri Abhyankar 78124121865SAdrian Maldonado /*@ 78224121865SAdrian Maldonado DMNetworkGetEdgeOffset - Get the offset for accessing the variable associated with the given edge from the local subvector. 78324121865SAdrian Maldonado 78424121865SAdrian Maldonado Not Collective 78524121865SAdrian Maldonado 78624121865SAdrian Maldonado Input Parameters: 78724121865SAdrian Maldonado + dm - The DMNetwork object 78824121865SAdrian Maldonado - p - the edge point 78924121865SAdrian Maldonado 79024121865SAdrian Maldonado Output Parameters: 79124121865SAdrian Maldonado . offset - the offset 79224121865SAdrian Maldonado 79324121865SAdrian Maldonado Level: intermediate 79424121865SAdrian Maldonado 79524121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 79624121865SAdrian Maldonado @*/ 79724121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset) 79824121865SAdrian Maldonado { 79924121865SAdrian Maldonado PetscErrorCode ierr; 80024121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 80124121865SAdrian Maldonado 80224121865SAdrian Maldonado PetscFunctionBegin; 80324121865SAdrian Maldonado 80424121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr); 80524121865SAdrian Maldonado PetscFunctionReturn(0); 80624121865SAdrian Maldonado } 80724121865SAdrian Maldonado 80824121865SAdrian Maldonado /*@ 80924121865SAdrian Maldonado DMNetworkGetVertexOffset - Get the offset for accessing the variable associated with the given vertex from the local subvector. 81024121865SAdrian Maldonado 81124121865SAdrian Maldonado Not Collective 81224121865SAdrian Maldonado 81324121865SAdrian Maldonado Input Parameters: 81424121865SAdrian Maldonado + dm - The DMNetwork object 81524121865SAdrian Maldonado - p - the vertex point 81624121865SAdrian Maldonado 81724121865SAdrian Maldonado Output Parameters: 81824121865SAdrian Maldonado . offset - the offset 81924121865SAdrian Maldonado 82024121865SAdrian Maldonado Level: intermediate 82124121865SAdrian Maldonado 82224121865SAdrian Maldonado .seealso: DMNetworkGetVariableGlobalOffset, DMGetLocalVector 82324121865SAdrian Maldonado @*/ 82424121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset) 82524121865SAdrian Maldonado { 82624121865SAdrian Maldonado PetscErrorCode ierr; 82724121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 82824121865SAdrian Maldonado 82924121865SAdrian Maldonado PetscFunctionBegin; 83024121865SAdrian Maldonado 83124121865SAdrian Maldonado p -= network->vStart; 83224121865SAdrian Maldonado 83324121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr); 83424121865SAdrian Maldonado PetscFunctionReturn(0); 83524121865SAdrian Maldonado } 8365f2c45f1SShri Abhyankar /*@ 8375f2c45f1SShri Abhyankar DMNetworkAddNumVariables - Add number of variables associated with a given point. 8385f2c45f1SShri Abhyankar 8395f2c45f1SShri Abhyankar Not Collective 8405f2c45f1SShri Abhyankar 8415f2c45f1SShri Abhyankar Input Parameters: 8425f2c45f1SShri Abhyankar + dm - The DMNetworkObject 8435f2c45f1SShri Abhyankar . p - the vertex/edge point 8445f2c45f1SShri Abhyankar - nvar - number of additional variables 8455f2c45f1SShri Abhyankar 8465f2c45f1SShri Abhyankar Level: intermediate 8475f2c45f1SShri Abhyankar 8485f2c45f1SShri Abhyankar .seealso: DMNetworkSetNumVariables 8495f2c45f1SShri Abhyankar @*/ 8505f2c45f1SShri Abhyankar PetscErrorCode DMNetworkAddNumVariables(DM dm,PetscInt p,PetscInt nvar) 8515f2c45f1SShri Abhyankar { 8525f2c45f1SShri Abhyankar PetscErrorCode ierr; 8535f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 8545f2c45f1SShri Abhyankar 8555f2c45f1SShri Abhyankar PetscFunctionBegin; 8565f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 8575f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8585f2c45f1SShri Abhyankar } 8595f2c45f1SShri Abhyankar 86027f51fceSHong Zhang /*@ 86127f51fceSHong Zhang DMNetworkGetNumVariables - Gets number of variables for a vertex/edge point. 86227f51fceSHong Zhang 86327f51fceSHong Zhang Not Collective 86427f51fceSHong Zhang 86527f51fceSHong Zhang Input Parameters: 86627f51fceSHong Zhang + dm - The DMNetworkObject 86727f51fceSHong Zhang - p - the vertex/edge point 86827f51fceSHong Zhang 86927f51fceSHong Zhang Output Parameters: 87027f51fceSHong Zhang . nvar - number of variables 87127f51fceSHong Zhang 87227f51fceSHong Zhang Level: intermediate 87327f51fceSHong Zhang 87427f51fceSHong Zhang .seealso: DMNetworkAddNumVariables, DMNetworkSddNumVariables 87527f51fceSHong Zhang @*/ 87627f51fceSHong Zhang PetscErrorCode DMNetworkGetNumVariables(DM dm,PetscInt p,PetscInt *nvar) 87727f51fceSHong Zhang { 87827f51fceSHong Zhang PetscErrorCode ierr; 87927f51fceSHong Zhang DM_Network *network = (DM_Network*)dm->data; 88027f51fceSHong Zhang 88127f51fceSHong Zhang PetscFunctionBegin; 88227f51fceSHong Zhang ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 88327f51fceSHong Zhang PetscFunctionReturn(0); 88427f51fceSHong Zhang } 88527f51fceSHong Zhang 8865f2c45f1SShri Abhyankar /*@ 8875f2c45f1SShri Abhyankar DMNetworkSetNumVariables - Sets number of variables for a vertex/edge point. 8885f2c45f1SShri Abhyankar 8895f2c45f1SShri Abhyankar Not Collective 8905f2c45f1SShri Abhyankar 8915f2c45f1SShri Abhyankar Input Parameters: 8925f2c45f1SShri Abhyankar + dm - The DMNetworkObject 8935f2c45f1SShri Abhyankar . p - the vertex/edge point 8945f2c45f1SShri Abhyankar - nvar - number of variables 8955f2c45f1SShri Abhyankar 8965f2c45f1SShri Abhyankar Level: intermediate 8975f2c45f1SShri Abhyankar 8985f2c45f1SShri Abhyankar .seealso: DMNetworkAddNumVariables 8995f2c45f1SShri Abhyankar @*/ 9005f2c45f1SShri Abhyankar PetscErrorCode DMNetworkSetNumVariables(DM dm,PetscInt p,PetscInt nvar) 9015f2c45f1SShri Abhyankar { 9025f2c45f1SShri Abhyankar PetscErrorCode ierr; 9035f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9045f2c45f1SShri Abhyankar 9055f2c45f1SShri Abhyankar PetscFunctionBegin; 9065f2c45f1SShri Abhyankar ierr = PetscSectionSetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 9075f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9085f2c45f1SShri Abhyankar } 9095f2c45f1SShri Abhyankar 9105f2c45f1SShri Abhyankar /* Sets up the array that holds the data for all components and its associated section. This 9115f2c45f1SShri Abhyankar function is called during DMSetUp() */ 9125f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 9135f2c45f1SShri Abhyankar { 9145f2c45f1SShri Abhyankar PetscErrorCode ierr; 9155f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 916e53b5ba3SHong Zhang PetscInt arr_size,p,offset,offsetp,ncomp,i; 9175f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 9185f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 9195f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 9205f2c45f1SShri Abhyankar 9215f2c45f1SShri Abhyankar PetscFunctionBegin; 9225f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 9235f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 92475b160a0SShri Abhyankar ierr = PetscMalloc1(arr_size,&network->componentdataarray);CHKERRQ(ierr); 9255f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 9265f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 9275f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 9285f2c45f1SShri Abhyankar /* Copy header */ 9295f2c45f1SShri Abhyankar header = &network->header[p]; 930302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offsetp,header,network->dataheadersize*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 9315f2c45f1SShri Abhyankar /* Copy data */ 9325f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 9335f2c45f1SShri Abhyankar ncomp = header->ndata; 9345f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 9355f2c45f1SShri Abhyankar offset = offsetp + network->dataheadersize + header->offset[i]; 936302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 9375f2c45f1SShri Abhyankar } 9385f2c45f1SShri Abhyankar } 9395f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9405f2c45f1SShri Abhyankar } 9415f2c45f1SShri Abhyankar 9425f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 9435f2c45f1SShri Abhyankar PetscErrorCode DMNetworkVariablesSetUp(DM dm) 9445f2c45f1SShri Abhyankar { 9455f2c45f1SShri Abhyankar PetscErrorCode ierr; 9465f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9475f2c45f1SShri Abhyankar 9485f2c45f1SShri Abhyankar PetscFunctionBegin; 9495f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 9505f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9515f2c45f1SShri Abhyankar } 9525f2c45f1SShri Abhyankar 9535f2c45f1SShri Abhyankar /*@C 9545f2c45f1SShri Abhyankar DMNetworkGetComponentDataArray - Returns the component data array 9555f2c45f1SShri Abhyankar 9565f2c45f1SShri Abhyankar Not Collective 9575f2c45f1SShri Abhyankar 9585f2c45f1SShri Abhyankar Input Parameters: 9595f2c45f1SShri Abhyankar . dm - The DMNetwork Object 9605f2c45f1SShri Abhyankar 9615f2c45f1SShri Abhyankar Output Parameters: 9625f2c45f1SShri Abhyankar . componentdataarray - array that holds data for all components 9635f2c45f1SShri Abhyankar 9645f2c45f1SShri Abhyankar Level: intermediate 9655f2c45f1SShri Abhyankar 966a730d845SHong Zhang .seealso: DMNetworkGetComponentKeyOffset, DMNetworkGetNumComponents 9675f2c45f1SShri Abhyankar @*/ 9685f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetComponentDataArray(DM dm,DMNetworkComponentGenericDataType **componentdataarray) 9695f2c45f1SShri Abhyankar { 9705f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9715f2c45f1SShri Abhyankar 9725f2c45f1SShri Abhyankar PetscFunctionBegin; 9735f2c45f1SShri Abhyankar *componentdataarray = network->componentdataarray; 9745f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9755f2c45f1SShri Abhyankar } 9765f2c45f1SShri Abhyankar 97724121865SAdrian Maldonado /* Get a subsection from a range of points */ 97824121865SAdrian Maldonado PetscErrorCode DMNetworkGetSubSection_private(PetscSection master, PetscInt pstart, PetscInt pend,PetscSection *subsection) 97924121865SAdrian Maldonado { 98024121865SAdrian Maldonado PetscErrorCode ierr; 98124121865SAdrian Maldonado PetscInt i, nvar; 98224121865SAdrian Maldonado 98324121865SAdrian Maldonado PetscFunctionBegin; 98424121865SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)master), subsection);CHKERRQ(ierr); 98524121865SAdrian Maldonado ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr); 98624121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 98724121865SAdrian Maldonado ierr = PetscSectionGetDof(master,i,&nvar);CHKERRQ(ierr); 98824121865SAdrian Maldonado ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr); 98924121865SAdrian Maldonado } 99024121865SAdrian Maldonado 99124121865SAdrian Maldonado ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr); 99224121865SAdrian Maldonado PetscFunctionReturn(0); 99324121865SAdrian Maldonado } 99424121865SAdrian Maldonado 99524121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */ 99624121865SAdrian Maldonado PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 99724121865SAdrian Maldonado { 99824121865SAdrian Maldonado PetscErrorCode ierr; 99924121865SAdrian Maldonado PetscInt i, *subpoints; 100024121865SAdrian Maldonado 100124121865SAdrian Maldonado PetscFunctionBegin; 100224121865SAdrian Maldonado /* Create index sets to map from "points" to "subpoints" */ 100324121865SAdrian Maldonado ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr); 100424121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 100524121865SAdrian Maldonado subpoints[i - pstart] = i; 100624121865SAdrian Maldonado } 1007459726d8SSatish Balay ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr); 100824121865SAdrian Maldonado ierr = PetscFree(subpoints);CHKERRQ(ierr); 100924121865SAdrian Maldonado PetscFunctionReturn(0); 101024121865SAdrian Maldonado } 101124121865SAdrian Maldonado 101224121865SAdrian Maldonado /*@ 101324121865SAdrian Maldonado DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute. 101424121865SAdrian Maldonado 101524121865SAdrian Maldonado Collective 101624121865SAdrian Maldonado 101724121865SAdrian Maldonado Input Parameters: 101824121865SAdrian Maldonado . dm - The DMNetworkObject 101924121865SAdrian Maldonado 102024121865SAdrian Maldonado Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 102124121865SAdrian Maldonado 102224121865SAdrian Maldonado points = [0 1 2 3 4 5 6] 102324121865SAdrian Maldonado 1024caf410d2SHong Zhang where edges = [0,1,2,3] and vertices = [4,5,6]. The new orderings will be specific to the subset (i.e vertices = [0,1,2] <- [4,5,6]). 102524121865SAdrian Maldonado 102624121865SAdrian Maldonado With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 102724121865SAdrian Maldonado 102824121865SAdrian Maldonado Level: intermediate 102924121865SAdrian Maldonado 103024121865SAdrian Maldonado @*/ 103124121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 103224121865SAdrian Maldonado { 103324121865SAdrian Maldonado PetscErrorCode ierr; 103424121865SAdrian Maldonado MPI_Comm comm; 10359852e123SBarry Smith PetscMPIInt rank, size; 103624121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 103724121865SAdrian Maldonado 1038eab1376dSHong Zhang PetscFunctionBegin; 103924121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 104024121865SAdrian Maldonado ierr = MPI_Comm_rank(comm, &rank);CHKERRQ(ierr); 10419852e123SBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 104224121865SAdrian Maldonado 104324121865SAdrian Maldonado /* Create maps for vertices and edges */ 104424121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 104524121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr); 104624121865SAdrian Maldonado 104724121865SAdrian Maldonado /* Create local sub-sections */ 104824121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr); 104924121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr); 105024121865SAdrian Maldonado 10519852e123SBarry Smith if (size > 1) { 105224121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 105324121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr); 105424121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr); 105524121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr); 105624121865SAdrian Maldonado } else { 105724121865SAdrian Maldonado /* create structures for vertex */ 105824121865SAdrian Maldonado ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr); 105924121865SAdrian Maldonado /* create structures for edge */ 106024121865SAdrian Maldonado ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr); 106124121865SAdrian Maldonado } 106224121865SAdrian Maldonado 106324121865SAdrian Maldonado 106424121865SAdrian Maldonado /* Add viewers */ 106524121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr); 106624121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr); 106724121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr); 106824121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr); 106924121865SAdrian Maldonado 107024121865SAdrian Maldonado PetscFunctionReturn(0); 107124121865SAdrian Maldonado } 10727b6afd5bSHong Zhang 10735f2c45f1SShri Abhyankar /*@ 10745f2c45f1SShri Abhyankar DMNetworkDistribute - Distributes the network and moves associated component data. 10755f2c45f1SShri Abhyankar 10765f2c45f1SShri Abhyankar Collective 10775f2c45f1SShri Abhyankar 10785f2c45f1SShri Abhyankar Input Parameter: 1079d3464fd4SAdrian Maldonado + DM - the DMNetwork object 10805f2c45f1SShri Abhyankar - overlap - The overlap of partitions, 0 is the default 10815f2c45f1SShri Abhyankar 10825f2c45f1SShri Abhyankar Notes: 10838b171c8eSHong Zhang Distributes the network with <overlap>-overlapping partitioning of the edges. 10845f2c45f1SShri Abhyankar 10855f2c45f1SShri Abhyankar Level: intermediate 10865f2c45f1SShri Abhyankar 10875f2c45f1SShri Abhyankar .seealso: DMNetworkCreate 10885f2c45f1SShri Abhyankar @*/ 1089d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap) 10905f2c45f1SShri Abhyankar { 1091d3464fd4SAdrian Maldonado MPI_Comm comm; 10925f2c45f1SShri Abhyankar PetscErrorCode ierr; 1093d3464fd4SAdrian Maldonado PetscMPIInt size; 1094d3464fd4SAdrian Maldonado DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data); 1095d3464fd4SAdrian Maldonado DM_Network *newDMnetwork; 1096caf410d2SHong Zhang PetscSF pointsf=NULL; 10975f2c45f1SShri Abhyankar DM newDM; 1098caf410d2SHong Zhang PetscInt j,e,v,offset,*subnetvtx; 109951ac5effSHong Zhang PetscPartitioner part; 1100b9c6e19dSShri Abhyankar DMNetworkComponentHeader header; 11015f2c45f1SShri Abhyankar 11025f2c45f1SShri Abhyankar PetscFunctionBegin; 1103d3464fd4SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr); 1104d3464fd4SAdrian Maldonado ierr = MPI_Comm_size(comm, &size);CHKERRQ(ierr); 1105d3464fd4SAdrian Maldonado if (size == 1) PetscFunctionReturn(0); 1106d3464fd4SAdrian Maldonado 1107d3464fd4SAdrian Maldonado ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr); 11085f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 11095f2c45f1SShri Abhyankar newDMnetwork->dataheadersize = sizeof(struct _p_DMNetworkComponentHeader)/sizeof(DMNetworkComponentGenericDataType); 111051ac5effSHong Zhang 111151ac5effSHong Zhang /* Enable runtime options for petscpartitioner */ 111251ac5effSHong Zhang ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr); 111351ac5effSHong Zhang ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 111451ac5effSHong Zhang 11155f2c45f1SShri Abhyankar /* Distribute plex dm and dof section */ 111680cf41d5SMatthew G. Knepley ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 111751ac5effSHong Zhang 11185f2c45f1SShri Abhyankar /* Distribute dof section */ 1119d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DofSection);CHKERRQ(ierr); 11205f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 1121d3464fd4SAdrian Maldonado ierr = PetscSectionCreate(PetscObjectComm((PetscObject)*dm),&newDMnetwork->DataSection);CHKERRQ(ierr); 112251ac5effSHong Zhang 11235f2c45f1SShri Abhyankar /* Distribute data and associated section */ 112431da1fc8SHong Zhang ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 112524121865SAdrian Maldonado 11265f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 11275f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 11285f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 11295f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 11306fefedf4SHong Zhang newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart; 11316fefedf4SHong Zhang newDMnetwork->NVertices = oldDMnetwork->NVertices; 11325f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 113324121865SAdrian Maldonado 11345f2c45f1SShri Abhyankar /* Set Dof section as the default section for dm */ 1135e87a4003SBarry Smith ierr = DMSetSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 1136e87a4003SBarry Smith ierr = DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 11375f2c45f1SShri Abhyankar 1138b9c6e19dSShri Abhyankar /* Set up subnetwork info in the newDM */ 1139b9c6e19dSShri Abhyankar newDMnetwork->nsubnet = oldDMnetwork->nsubnet; 1140caf410d2SHong Zhang newDMnetwork->ncsubnet = oldDMnetwork->ncsubnet; 1141b9c6e19dSShri Abhyankar ierr = PetscCalloc1(newDMnetwork->nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr); 1142b9c6e19dSShri Abhyankar /* Copy over the global number of vertices and edges in each subnetwork. Note that these are already 1143b9c6e19dSShri Abhyankar calculated in DMNetworkLayoutSetUp() 1144b9c6e19dSShri Abhyankar */ 1145b9c6e19dSShri Abhyankar for(j=0; j < newDMnetwork->nsubnet; j++) { 1146b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx; 1147b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge; 1148b9c6e19dSShri Abhyankar } 1149b9c6e19dSShri Abhyankar 1150b9c6e19dSShri Abhyankar for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) { 1151b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1152b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1153b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].nedge++; 1154b9c6e19dSShri Abhyankar } 1155b9c6e19dSShri Abhyankar 1156b9c6e19dSShri Abhyankar for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) { 1157b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1158b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1159b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].nvtx++; 1160b9c6e19dSShri Abhyankar } 1161b9c6e19dSShri Abhyankar 1162b9c6e19dSShri Abhyankar /* Now create the vertices and edge arrays for the subnetworks */ 1163caf410d2SHong Zhang ierr = PetscCalloc1(newDMnetwork->vEnd-newDMnetwork->vStart,&newDMnetwork->subnetvtx);CHKERRQ(ierr); 1164caf410d2SHong Zhang subnetvtx = newDMnetwork->subnetvtx; 1165caf410d2SHong Zhang 1166b9c6e19dSShri Abhyankar for (j=0; j<newDMnetwork->nsubnet; j++) { 1167b9c6e19dSShri Abhyankar ierr = PetscCalloc1(newDMnetwork->subnet[j].nedge,&newDMnetwork->subnet[j].edges);CHKERRQ(ierr); 1168caf410d2SHong Zhang newDMnetwork->subnet[j].vertices = subnetvtx; 1169caf410d2SHong Zhang subnetvtx += newDMnetwork->subnet[j].nvtx; 1170caf410d2SHong Zhang 1171b9c6e19dSShri Abhyankar /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. 1172b9c6e19dSShri Abhyankar These get updated when the vertices and edges are added. */ 1173b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0; 1174b9c6e19dSShri Abhyankar } 1175b9c6e19dSShri Abhyankar 1176b9c6e19dSShri Abhyankar /* Set the vertices and edges in each subnetwork */ 1177b9c6e19dSShri Abhyankar for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++ ) { 1178b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1179b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1180b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e; 1181b9c6e19dSShri Abhyankar } 1182b9c6e19dSShri Abhyankar 1183b9c6e19dSShri Abhyankar for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++ ) { 1184b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1185b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1186b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v; 1187b9c6e19dSShri Abhyankar } 1188b9c6e19dSShri Abhyankar 1189caf410d2SHong Zhang newDM->setupcalled = (*dm)->setupcalled; 1190caf410d2SHong Zhang 119124121865SAdrian Maldonado /* Destroy point SF */ 119224121865SAdrian Maldonado ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 119324121865SAdrian Maldonado 1194d3464fd4SAdrian Maldonado ierr = DMDestroy(dm);CHKERRQ(ierr); 1195d3464fd4SAdrian Maldonado *dm = newDM; 11965f2c45f1SShri Abhyankar PetscFunctionReturn(0); 11975f2c45f1SShri Abhyankar } 11985f2c45f1SShri Abhyankar 119924121865SAdrian Maldonado /*@C 120024121865SAdrian Maldonado PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering. 120124121865SAdrian Maldonado 120224121865SAdrian Maldonado Input Parameters: 120324121865SAdrian Maldonado + masterSF - the original SF structure 120424121865SAdrian Maldonado - map - a ISLocalToGlobal mapping that contains the subset of points 120524121865SAdrian Maldonado 120624121865SAdrian Maldonado Output Parameters: 120724121865SAdrian Maldonado . subSF - a subset of the masterSF for the desired subset. 120824121865SAdrian Maldonado */ 120924121865SAdrian Maldonado 121024121865SAdrian Maldonado PetscErrorCode PetscSFGetSubSF(PetscSF mastersf, ISLocalToGlobalMapping map, PetscSF *subSF) { 121124121865SAdrian Maldonado 121224121865SAdrian Maldonado PetscErrorCode ierr; 121324121865SAdrian Maldonado PetscInt nroots, nleaves, *ilocal_sub; 121424121865SAdrian Maldonado PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 121524121865SAdrian Maldonado PetscInt *local_points, *remote_points; 121624121865SAdrian Maldonado PetscSFNode *iremote_sub; 121724121865SAdrian Maldonado const PetscInt *ilocal; 121824121865SAdrian Maldonado const PetscSFNode *iremote; 121924121865SAdrian Maldonado 122024121865SAdrian Maldonado PetscFunctionBegin; 122124121865SAdrian Maldonado ierr = PetscSFGetGraph(mastersf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 122224121865SAdrian Maldonado 122324121865SAdrian Maldonado /* Look for leaves that pertain to the subset of points. Get the local ordering */ 122424121865SAdrian Maldonado ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr); 122524121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr); 122624121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 122724121865SAdrian Maldonado if (ilocal_map[i] != -1) nleaves_sub += 1; 122824121865SAdrian Maldonado } 122924121865SAdrian Maldonado /* Re-number ilocal with subset numbering. Need information from roots */ 123024121865SAdrian Maldonado ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr); 123124121865SAdrian Maldonado for (i = 0; i < nroots; i++) local_points[i] = i; 123224121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr); 123324121865SAdrian Maldonado ierr = PetscSFBcastBegin(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 123424121865SAdrian Maldonado ierr = PetscSFBcastEnd(mastersf, MPIU_INT, local_points, remote_points);CHKERRQ(ierr); 123524121865SAdrian Maldonado /* Fill up graph using local (that is, local to the subset) numbering. */ 12364b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr); 12374b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr); 123824121865SAdrian Maldonado nleaves_sub = 0; 123924121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 124024121865SAdrian Maldonado if (ilocal_map[i] != -1) { 124124121865SAdrian Maldonado ilocal_sub[nleaves_sub] = ilocal_map[i]; 12424b70a8deSAdrian Maldonado iremote_sub[nleaves_sub].rank = iremote[i].rank; 124324121865SAdrian Maldonado iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 124424121865SAdrian Maldonado nleaves_sub += 1; 124524121865SAdrian Maldonado } 124624121865SAdrian Maldonado } 124724121865SAdrian Maldonado ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr); 124824121865SAdrian Maldonado ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr); 124924121865SAdrian Maldonado 125024121865SAdrian Maldonado /* Create new subSF */ 125124121865SAdrian Maldonado ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr); 125224121865SAdrian Maldonado ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr); 12534b70a8deSAdrian Maldonado ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr); 125424121865SAdrian Maldonado ierr = PetscFree(ilocal_map);CHKERRQ(ierr); 12554b70a8deSAdrian Maldonado ierr = PetscFree(iremote_sub);CHKERRQ(ierr); 125624121865SAdrian Maldonado PetscFunctionReturn(0); 125724121865SAdrian Maldonado } 125824121865SAdrian Maldonado 12595f2c45f1SShri Abhyankar /*@C 12605f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 12615f2c45f1SShri Abhyankar 12625f2c45f1SShri Abhyankar Not Collective 12635f2c45f1SShri Abhyankar 12645f2c45f1SShri Abhyankar Input Parameters: 12655f2c45f1SShri Abhyankar + dm - The DMNetwork object 12665f2c45f1SShri Abhyankar - p - the vertex point 12675f2c45f1SShri Abhyankar 12685f2c45f1SShri Abhyankar Output Paramters: 12695f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 12705f2c45f1SShri Abhyankar - edges - List of edge points 12715f2c45f1SShri Abhyankar 12725f2c45f1SShri Abhyankar Level: intermediate 12735f2c45f1SShri Abhyankar 12745f2c45f1SShri Abhyankar Fortran Notes: 12755f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 12765f2c45f1SShri Abhyankar include petsc.h90 in your code. 12775f2c45f1SShri Abhyankar 1278d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices 12795f2c45f1SShri Abhyankar @*/ 12805f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 12815f2c45f1SShri Abhyankar { 12825f2c45f1SShri Abhyankar PetscErrorCode ierr; 12835f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 12845f2c45f1SShri Abhyankar 12855f2c45f1SShri Abhyankar PetscFunctionBegin; 12865f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 12875f2c45f1SShri Abhyankar ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr); 12885f2c45f1SShri Abhyankar PetscFunctionReturn(0); 12895f2c45f1SShri Abhyankar } 12905f2c45f1SShri Abhyankar 12915f2c45f1SShri Abhyankar /*@C 1292d842c372SHong Zhang DMNetworkGetConnectedVertices - Return the connected vertices for this edge point 12935f2c45f1SShri Abhyankar 12945f2c45f1SShri Abhyankar Not Collective 12955f2c45f1SShri Abhyankar 12965f2c45f1SShri Abhyankar Input Parameters: 12975f2c45f1SShri Abhyankar + dm - The DMNetwork object 12985f2c45f1SShri Abhyankar - p - the edge point 12995f2c45f1SShri Abhyankar 13005f2c45f1SShri Abhyankar Output Paramters: 13015f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 13025f2c45f1SShri Abhyankar 13035f2c45f1SShri Abhyankar Level: intermediate 13045f2c45f1SShri Abhyankar 13055f2c45f1SShri Abhyankar Fortran Notes: 13065f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 13075f2c45f1SShri Abhyankar include petsc.h90 in your code. 13085f2c45f1SShri Abhyankar 13095f2c45f1SShri Abhyankar .seealso: DMNetworkCreate, DMNetworkGetSupportingEdges 13105f2c45f1SShri Abhyankar @*/ 1311d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[]) 13125f2c45f1SShri Abhyankar { 13135f2c45f1SShri Abhyankar PetscErrorCode ierr; 13145f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 13155f2c45f1SShri Abhyankar 13165f2c45f1SShri Abhyankar PetscFunctionBegin; 13175f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 13185f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13195f2c45f1SShri Abhyankar } 13205f2c45f1SShri Abhyankar 13215f2c45f1SShri Abhyankar /*@ 13225f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 13235f2c45f1SShri Abhyankar 13245f2c45f1SShri Abhyankar Not Collective 13255f2c45f1SShri Abhyankar 13265f2c45f1SShri Abhyankar Input Parameters: 13275f2c45f1SShri Abhyankar + dm - The DMNetwork object 13285f2c45f1SShri Abhyankar . p - the vertex point 13295f2c45f1SShri Abhyankar 13305f2c45f1SShri Abhyankar Output Parameter: 13315f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 13325f2c45f1SShri Abhyankar 13335f2c45f1SShri Abhyankar Level: intermediate 13345f2c45f1SShri Abhyankar 1335d842c372SHong Zhang .seealso: DMNetworkCreate, DMNetworkGetConnectedVertices, DMNetworkGetVertexRange 13365f2c45f1SShri Abhyankar @*/ 13375f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 13385f2c45f1SShri Abhyankar { 13395f2c45f1SShri Abhyankar PetscErrorCode ierr; 13405f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 13415f2c45f1SShri Abhyankar PetscInt offsetg; 13425f2c45f1SShri Abhyankar PetscSection sectiong; 13435f2c45f1SShri Abhyankar 13445f2c45f1SShri Abhyankar PetscFunctionBegin; 1345caf410d2SHong Zhang if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first"); 13465f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 1347e87a4003SBarry Smith ierr = DMGetGlobalSection(network->plex,§iong);CHKERRQ(ierr); 13485f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 13495f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 13505f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13515f2c45f1SShri Abhyankar } 13525f2c45f1SShri Abhyankar 13535f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 13545f2c45f1SShri Abhyankar { 13555f2c45f1SShri Abhyankar PetscErrorCode ierr; 13565f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 13575f2c45f1SShri Abhyankar 13585f2c45f1SShri Abhyankar PetscFunctionBegin; 13595f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 13605f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 13615f2c45f1SShri Abhyankar 1362e87a4003SBarry Smith ierr = DMSetSection(network->plex,network->DofSection);CHKERRQ(ierr); 1363e87a4003SBarry Smith ierr = DMGetGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 1364*9e1d080bSHong Zhang 1365*9e1d080bSHong Zhang dm->setupcalled = PETSC_TRUE; 1366*9e1d080bSHong Zhang ierr = DMViewFromOptions(dm,NULL,"-dm_view");CHKERRQ(ierr); 13675f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13685f2c45f1SShri Abhyankar } 13695f2c45f1SShri Abhyankar 13701ad426b7SHong Zhang /*@ 137117df6e9eSHong Zhang DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 13721ad426b7SHong Zhang -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 13731ad426b7SHong Zhang 13741ad426b7SHong Zhang Collective 13751ad426b7SHong Zhang 13761ad426b7SHong Zhang Input Parameters: 137783b2e829SHong Zhang + dm - The DMNetwork object 137883b2e829SHong Zhang . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 137983b2e829SHong Zhang - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 13801ad426b7SHong Zhang 13811ad426b7SHong Zhang Level: intermediate 13821ad426b7SHong Zhang 13831ad426b7SHong Zhang @*/ 138483b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 13851ad426b7SHong Zhang { 13861ad426b7SHong Zhang DM_Network *network=(DM_Network*)dm->data; 13878675203cSHong Zhang PetscErrorCode ierr; 138866b4e0ffSHong Zhang PetscInt nVertices = network->nVertices; 13891ad426b7SHong Zhang 13901ad426b7SHong Zhang PetscFunctionBegin; 139183b2e829SHong Zhang network->userEdgeJacobian = eflg; 139283b2e829SHong Zhang network->userVertexJacobian = vflg; 13938675203cSHong Zhang 13948675203cSHong Zhang if (eflg && !network->Je) { 13958675203cSHong Zhang ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 13968675203cSHong Zhang } 13978675203cSHong Zhang 139866b4e0ffSHong Zhang if (vflg && !network->Jv && nVertices) { 13998675203cSHong Zhang PetscInt i,*vptr,nedges,vStart=network->vStart; 140066b4e0ffSHong Zhang PetscInt nedges_total; 14018675203cSHong Zhang const PetscInt *edges; 14028675203cSHong Zhang 14038675203cSHong Zhang /* count nvertex_total */ 14048675203cSHong Zhang nedges_total = 0; 14058675203cSHong Zhang ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr); 14068675203cSHong Zhang 14078675203cSHong Zhang vptr[0] = 0; 14088675203cSHong Zhang for (i=0; i<nVertices; i++) { 14098675203cSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 14108675203cSHong Zhang nedges_total += nedges; 14118675203cSHong Zhang vptr[i+1] = vptr[i] + 2*nedges + 1; 14128675203cSHong Zhang } 14138675203cSHong Zhang 14148675203cSHong Zhang ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr); 14158675203cSHong Zhang network->Jvptr = vptr; 14168675203cSHong Zhang } 14171ad426b7SHong Zhang PetscFunctionReturn(0); 14181ad426b7SHong Zhang } 14191ad426b7SHong Zhang 14201ad426b7SHong Zhang /*@ 142183b2e829SHong Zhang DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 142283b2e829SHong Zhang 142383b2e829SHong Zhang Not Collective 142483b2e829SHong Zhang 142583b2e829SHong Zhang Input Parameters: 142683b2e829SHong Zhang + dm - The DMNetwork object 142783b2e829SHong Zhang . p - the edge point 14283e97b6e8SHong Zhang - J - array (size = 3) of Jacobian submatrices for this edge point: 14293e97b6e8SHong Zhang J[0]: this edge 1430d842c372SHong Zhang J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices() 143183b2e829SHong Zhang 143283b2e829SHong Zhang Level: intermediate 143383b2e829SHong Zhang 143483b2e829SHong Zhang .seealso: DMNetworkVertexSetMatrix 143583b2e829SHong Zhang @*/ 143683b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 143783b2e829SHong Zhang { 143883b2e829SHong Zhang DM_Network *network=(DM_Network*)dm->data; 143983b2e829SHong Zhang 144083b2e829SHong Zhang PetscFunctionBegin; 14418675203cSHong Zhang if (!network->Je) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 14428675203cSHong Zhang 14438675203cSHong Zhang if (J) { 1444883e35e8SHong Zhang network->Je[3*p] = J[0]; 1445883e35e8SHong Zhang network->Je[3*p+1] = J[1]; 1446883e35e8SHong Zhang network->Je[3*p+2] = J[2]; 14478675203cSHong Zhang } 144883b2e829SHong Zhang PetscFunctionReturn(0); 144983b2e829SHong Zhang } 145083b2e829SHong Zhang 145183b2e829SHong Zhang /*@ 145276ddfea5SHong Zhang DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 14531ad426b7SHong Zhang 14541ad426b7SHong Zhang Not Collective 14551ad426b7SHong Zhang 14561ad426b7SHong Zhang Input Parameters: 14571ad426b7SHong Zhang + dm - The DMNetwork object 14581ad426b7SHong Zhang . p - the vertex point 14593e97b6e8SHong Zhang - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 14603e97b6e8SHong Zhang J[0]: this vertex 14613e97b6e8SHong Zhang J[1+2*i]: i-th supporting edge 14623e97b6e8SHong Zhang J[1+2*i+1]: i-th connected vertex 14631ad426b7SHong Zhang 14641ad426b7SHong Zhang Level: intermediate 14651ad426b7SHong Zhang 146683b2e829SHong Zhang .seealso: DMNetworkEdgeSetMatrix 14671ad426b7SHong Zhang @*/ 1468883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 14695f2c45f1SShri Abhyankar { 14705f2c45f1SShri Abhyankar PetscErrorCode ierr; 14715f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 14728675203cSHong Zhang PetscInt i,*vptr,nedges,vStart=network->vStart; 1473883e35e8SHong Zhang const PetscInt *edges; 14745f2c45f1SShri Abhyankar 14755f2c45f1SShri Abhyankar PetscFunctionBegin; 14768675203cSHong Zhang if (!network->Jv) SETERRQ(PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 1477883e35e8SHong Zhang 14788675203cSHong Zhang if (J) { 1479883e35e8SHong Zhang vptr = network->Jvptr; 14803e97b6e8SHong Zhang network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 14813e97b6e8SHong Zhang 14823e97b6e8SHong Zhang /* Set Jacobian for each supporting edge and connected vertex */ 1483883e35e8SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 1484883e35e8SHong Zhang for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 14858675203cSHong Zhang } 1486883e35e8SHong Zhang PetscFunctionReturn(0); 1487883e35e8SHong Zhang } 1488883e35e8SHong Zhang 1489e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 14905cf7da58SHong Zhang { 14915cf7da58SHong Zhang PetscErrorCode ierr; 14925cf7da58SHong Zhang PetscInt j; 14935cf7da58SHong Zhang PetscScalar val=(PetscScalar)ncols; 14945cf7da58SHong Zhang 14955cf7da58SHong Zhang PetscFunctionBegin; 14965cf7da58SHong Zhang if (!ghost) { 14975cf7da58SHong Zhang for (j=0; j<nrows; j++) { 14985cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 14995cf7da58SHong Zhang } 15005cf7da58SHong Zhang } else { 15015cf7da58SHong Zhang for (j=0; j<nrows; j++) { 15025cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 15035cf7da58SHong Zhang } 15045cf7da58SHong Zhang } 15055cf7da58SHong Zhang PetscFunctionReturn(0); 15065cf7da58SHong Zhang } 15075cf7da58SHong Zhang 1508e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 15095cf7da58SHong Zhang { 15105cf7da58SHong Zhang PetscErrorCode ierr; 15115cf7da58SHong Zhang PetscInt j,ncols_u; 15125cf7da58SHong Zhang PetscScalar val; 15135cf7da58SHong Zhang 15145cf7da58SHong Zhang PetscFunctionBegin; 15155cf7da58SHong Zhang if (!ghost) { 15165cf7da58SHong Zhang for (j=0; j<nrows; j++) { 15175cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 15185cf7da58SHong Zhang val = (PetscScalar)ncols_u; 15195cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 15205cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 15215cf7da58SHong Zhang } 15225cf7da58SHong Zhang } else { 15235cf7da58SHong Zhang for (j=0; j<nrows; j++) { 15245cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 15255cf7da58SHong Zhang val = (PetscScalar)ncols_u; 15265cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 15275cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 15285cf7da58SHong Zhang } 15295cf7da58SHong Zhang } 15305cf7da58SHong Zhang PetscFunctionReturn(0); 15315cf7da58SHong Zhang } 15325cf7da58SHong Zhang 1533e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 15345cf7da58SHong Zhang { 15355cf7da58SHong Zhang PetscErrorCode ierr; 15365cf7da58SHong Zhang 15375cf7da58SHong Zhang PetscFunctionBegin; 15385cf7da58SHong Zhang if (Ju) { 15395cf7da58SHong Zhang ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 15405cf7da58SHong Zhang } else { 15415cf7da58SHong Zhang ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 15425cf7da58SHong Zhang } 15435cf7da58SHong Zhang PetscFunctionReturn(0); 15445cf7da58SHong Zhang } 15455cf7da58SHong Zhang 1546e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1547883e35e8SHong Zhang { 1548883e35e8SHong Zhang PetscErrorCode ierr; 1549883e35e8SHong Zhang PetscInt j,*cols; 1550883e35e8SHong Zhang PetscScalar *zeros; 1551883e35e8SHong Zhang 1552883e35e8SHong Zhang PetscFunctionBegin; 1553883e35e8SHong Zhang ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 1554883e35e8SHong Zhang for (j=0; j<ncols; j++) cols[j] = j+ cstart; 1555883e35e8SHong Zhang ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 1556883e35e8SHong Zhang ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 15571ad426b7SHong Zhang PetscFunctionReturn(0); 15581ad426b7SHong Zhang } 1559a4e85ca8SHong Zhang 1560e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 15613e97b6e8SHong Zhang { 15623e97b6e8SHong Zhang PetscErrorCode ierr; 15633e97b6e8SHong Zhang PetscInt j,M,N,row,col,ncols_u; 15643e97b6e8SHong Zhang const PetscInt *cols; 15653e97b6e8SHong Zhang PetscScalar zero=0.0; 15663e97b6e8SHong Zhang 15673e97b6e8SHong Zhang PetscFunctionBegin; 15683e97b6e8SHong Zhang ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 15693e97b6e8SHong 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); 15703e97b6e8SHong Zhang 15713e97b6e8SHong Zhang for (row=0; row<nrows; row++) { 15723e97b6e8SHong Zhang ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 15733e97b6e8SHong Zhang for (j=0; j<ncols_u; j++) { 15743e97b6e8SHong Zhang col = cols[j] + cstart; 15753e97b6e8SHong Zhang ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 15763e97b6e8SHong Zhang } 15773e97b6e8SHong Zhang ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 15783e97b6e8SHong Zhang } 15793e97b6e8SHong Zhang PetscFunctionReturn(0); 15803e97b6e8SHong Zhang } 15811ad426b7SHong Zhang 1582e0f69777SHong Zhang PETSC_STATIC_INLINE PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 1583a4e85ca8SHong Zhang { 1584a4e85ca8SHong Zhang PetscErrorCode ierr; 1585f4431b8cSHong Zhang 1586a4e85ca8SHong Zhang PetscFunctionBegin; 1587a4e85ca8SHong Zhang if (Ju) { 1588a4e85ca8SHong Zhang ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1589a4e85ca8SHong Zhang } else { 1590a4e85ca8SHong Zhang ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1591a4e85ca8SHong Zhang } 1592a4e85ca8SHong Zhang PetscFunctionReturn(0); 1593a4e85ca8SHong Zhang } 1594a4e85ca8SHong Zhang 159524121865SAdrian 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. 159624121865SAdrian Maldonado */ 159724121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 159824121865SAdrian Maldonado { 159924121865SAdrian Maldonado PetscErrorCode ierr; 160024121865SAdrian Maldonado PetscInt i,size,dof; 160124121865SAdrian Maldonado PetscInt *glob2loc; 160224121865SAdrian Maldonado 160324121865SAdrian Maldonado PetscFunctionBegin; 160424121865SAdrian Maldonado ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr); 160524121865SAdrian Maldonado ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr); 160624121865SAdrian Maldonado 160724121865SAdrian Maldonado for (i = 0; i < size; i++) { 160824121865SAdrian Maldonado ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr); 160924121865SAdrian Maldonado dof = (dof >= 0) ? dof : -(dof + 1); 161024121865SAdrian Maldonado glob2loc[i] = dof; 161124121865SAdrian Maldonado } 161224121865SAdrian Maldonado 161324121865SAdrian Maldonado ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 161424121865SAdrian Maldonado #if 0 161524121865SAdrian Maldonado ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 161624121865SAdrian Maldonado #endif 161724121865SAdrian Maldonado PetscFunctionReturn(0); 161824121865SAdrian Maldonado } 161924121865SAdrian Maldonado 162001ad2aeeSHong Zhang #include <petsc/private/matimpl.h> 1621*9e1d080bSHong Zhang 1622*9e1d080bSHong Zhang PetscErrorCode DMCreateMatrix_Network_Nest(DM dm,Mat *J) 16231ad426b7SHong Zhang { 16241ad426b7SHong Zhang PetscErrorCode ierr; 16251ad426b7SHong Zhang DM_Network *network = (DM_Network*)dm->data; 1626*9e1d080bSHong Zhang PetscMPIInt rank, size; 162724121865SAdrian Maldonado PetscInt eDof,vDof; 162824121865SAdrian Maldonado Mat j11,j12,j21,j22,bA[2][2]; 1629*9e1d080bSHong Zhang MPI_Comm comm; 163024121865SAdrian Maldonado ISLocalToGlobalMapping eISMap,vISMap; 163124121865SAdrian Maldonado 1632*9e1d080bSHong Zhang PetscFunctionBegin; 163324121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 163424121865SAdrian Maldonado ierr = MPI_Comm_rank(comm,&rank);CHKERRQ(ierr); 163524121865SAdrian Maldonado ierr = MPI_Comm_size(comm,&size);CHKERRQ(ierr); 163624121865SAdrian Maldonado 163724121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr); 163824121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr); 163924121865SAdrian Maldonado 164001ad2aeeSHong Zhang ierr = MatCreate(comm, &j11);CHKERRQ(ierr); 164124121865SAdrian Maldonado ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 164224121865SAdrian Maldonado ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr); 164324121865SAdrian Maldonado 164401ad2aeeSHong Zhang ierr = MatCreate(comm, &j12);CHKERRQ(ierr); 164524121865SAdrian Maldonado ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr); 164624121865SAdrian Maldonado ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr); 164724121865SAdrian Maldonado 164801ad2aeeSHong Zhang ierr = MatCreate(comm, &j21);CHKERRQ(ierr); 164924121865SAdrian Maldonado ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 165024121865SAdrian Maldonado ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr); 165124121865SAdrian Maldonado 165201ad2aeeSHong Zhang ierr = MatCreate(comm, &j22);CHKERRQ(ierr); 165324121865SAdrian Maldonado ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 165424121865SAdrian Maldonado ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr); 165524121865SAdrian Maldonado 16563f6a6bdaSHong Zhang bA[0][0] = j11; 16573f6a6bdaSHong Zhang bA[0][1] = j12; 16583f6a6bdaSHong Zhang bA[1][0] = j21; 16593f6a6bdaSHong Zhang bA[1][1] = j22; 166024121865SAdrian Maldonado 166124121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr); 166224121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr); 166324121865SAdrian Maldonado 166424121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr); 166524121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr); 166624121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr); 166724121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr); 166824121865SAdrian Maldonado 166924121865SAdrian Maldonado ierr = MatSetUp(j11);CHKERRQ(ierr); 167024121865SAdrian Maldonado ierr = MatSetUp(j12);CHKERRQ(ierr); 167124121865SAdrian Maldonado ierr = MatSetUp(j21);CHKERRQ(ierr); 167224121865SAdrian Maldonado ierr = MatSetUp(j22);CHKERRQ(ierr); 167324121865SAdrian Maldonado 167401ad2aeeSHong Zhang ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr); 167524121865SAdrian Maldonado ierr = MatSetUp(*J);CHKERRQ(ierr); 167624121865SAdrian Maldonado ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr); 167724121865SAdrian Maldonado ierr = MatDestroy(&j11);CHKERRQ(ierr); 167824121865SAdrian Maldonado ierr = MatDestroy(&j12);CHKERRQ(ierr); 167924121865SAdrian Maldonado ierr = MatDestroy(&j21);CHKERRQ(ierr); 168024121865SAdrian Maldonado ierr = MatDestroy(&j22);CHKERRQ(ierr); 168124121865SAdrian Maldonado 168224121865SAdrian Maldonado ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 168324121865SAdrian Maldonado ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1684*9e1d080bSHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 168524121865SAdrian Maldonado 168624121865SAdrian Maldonado /* Free structures */ 168724121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr); 168824121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr); 168924121865SAdrian Maldonado PetscFunctionReturn(0); 1690*9e1d080bSHong Zhang } 1691*9e1d080bSHong Zhang 1692*9e1d080bSHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 1693*9e1d080bSHong Zhang { 1694*9e1d080bSHong Zhang PetscErrorCode ierr; 1695*9e1d080bSHong Zhang DM_Network *network = (DM_Network*)dm->data; 1696*9e1d080bSHong Zhang PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 1697*9e1d080bSHong Zhang PetscInt cstart,ncols,j,e,v; 1698*9e1d080bSHong Zhang PetscBool ghost,ghost_vc,ghost2,isNest; 1699*9e1d080bSHong Zhang Mat Juser; 1700*9e1d080bSHong Zhang PetscSection sectionGlobal; 1701*9e1d080bSHong Zhang PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 1702*9e1d080bSHong Zhang const PetscInt *edges,*cone; 1703*9e1d080bSHong Zhang MPI_Comm comm; 1704*9e1d080bSHong Zhang MatType mtype; 1705*9e1d080bSHong Zhang Vec vd_nz,vo_nz; 1706*9e1d080bSHong Zhang PetscInt *dnnz,*onnz; 1707*9e1d080bSHong Zhang PetscScalar *vdnz,*vonz; 1708*9e1d080bSHong Zhang 1709*9e1d080bSHong Zhang PetscFunctionBegin; 1710*9e1d080bSHong Zhang mtype = dm->mattype; 1711*9e1d080bSHong Zhang ierr = PetscStrcmp(mtype,MATNEST,&isNest);CHKERRQ(ierr); 1712*9e1d080bSHong Zhang if (isNest) { 1713*9e1d080bSHong Zhang ierr = DMCreateMatrix_Network_Nest(dm,J);CHKERRQ(ierr); 1714*9e1d080bSHong Zhang PetscFunctionReturn(0); 1715*9e1d080bSHong Zhang } 1716*9e1d080bSHong Zhang 1717*9e1d080bSHong Zhang if (!network->userEdgeJacobian && !network->userVertexJacobian) { 1718a4e85ca8SHong Zhang /* user does not provide Jacobian blocks */ 1719*9e1d080bSHong Zhang ierr = DMCreateMatrix_Plex(network->plex,J);CHKERRQ(ierr); 1720*9e1d080bSHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 1721bfbc38dcSHong Zhang ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 17221ad426b7SHong Zhang PetscFunctionReturn(0); 17231ad426b7SHong Zhang } 17241ad426b7SHong Zhang 1725bfbc38dcSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 1726e87a4003SBarry Smith ierr = DMGetGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 1727bfbc38dcSHong Zhang ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 1728bfbc38dcSHong Zhang ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 17292a945128SHong Zhang 17302a945128SHong Zhang ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 17312a945128SHong Zhang ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 173289898e50SHong Zhang 173389898e50SHong Zhang /* (1) Set matrix preallocation */ 173489898e50SHong Zhang /*------------------------------*/ 1735840c2264SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1736840c2264SHong Zhang ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr); 1737840c2264SHong Zhang ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr); 1738840c2264SHong Zhang ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr); 1739840c2264SHong Zhang ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr); 1740840c2264SHong Zhang ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr); 1741840c2264SHong Zhang 174289898e50SHong Zhang /* Set preallocation for edges */ 174389898e50SHong Zhang /*-----------------------------*/ 1744840c2264SHong Zhang ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 1745840c2264SHong Zhang 1746bdcb62a2SHong Zhang ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr); 1747840c2264SHong Zhang for (e=eStart; e<eEnd; e++) { 1748840c2264SHong Zhang /* Get row indices */ 1749840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 1750840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 1751840c2264SHong Zhang if (nrows) { 1752840c2264SHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 1753840c2264SHong Zhang 17545cf7da58SHong Zhang /* Set preallocation for conntected vertices */ 1755d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1756840c2264SHong Zhang for (v=0; v<2; v++) { 1757840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 1758840c2264SHong Zhang 17598675203cSHong Zhang if (network->Je) { 1760840c2264SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 17618675203cSHong Zhang } else Juser = NULL; 1762840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 17635cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1764840c2264SHong Zhang } 1765840c2264SHong Zhang 176689898e50SHong Zhang /* Set preallocation for edge self */ 1767840c2264SHong Zhang cstart = rstart; 17688675203cSHong Zhang if (network->Je) { 1769840c2264SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 17708675203cSHong Zhang } else Juser = NULL; 17715cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1772840c2264SHong Zhang } 1773840c2264SHong Zhang } 1774840c2264SHong Zhang 177589898e50SHong Zhang /* Set preallocation for vertices */ 177689898e50SHong Zhang /*--------------------------------*/ 1777840c2264SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 17788675203cSHong Zhang if (vEnd - vStart) vptr = network->Jvptr; 1779840c2264SHong Zhang 1780840c2264SHong Zhang for (v=vStart; v<vEnd; v++) { 1781840c2264SHong Zhang /* Get row indices */ 1782840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1783840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 1784840c2264SHong Zhang if (!nrows) continue; 1785840c2264SHong Zhang 1786bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1787bdcb62a2SHong Zhang if (ghost) { 1788bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1789bdcb62a2SHong Zhang } else { 1790bdcb62a2SHong Zhang rows_v = rows; 1791bdcb62a2SHong Zhang } 1792bdcb62a2SHong Zhang 1793bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1794840c2264SHong Zhang 1795840c2264SHong Zhang /* Get supporting edges and connected vertices */ 1796840c2264SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1797840c2264SHong Zhang 1798840c2264SHong Zhang for (e=0; e<nedges; e++) { 1799840c2264SHong Zhang /* Supporting edges */ 1800840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1801840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1802840c2264SHong Zhang 18038675203cSHong Zhang if (network->Jv) { 1804840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 18058675203cSHong Zhang } else Juser = NULL; 1806bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 1807840c2264SHong Zhang 1808840c2264SHong Zhang /* Connected vertices */ 1809d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 1810840c2264SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 1811840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr); 1812840c2264SHong Zhang 1813840c2264SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1814840c2264SHong Zhang 18158675203cSHong Zhang if (network->Jv) { 1816840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 18178675203cSHong Zhang } else Juser = NULL; 1818e102a522SHong Zhang if (ghost_vc||ghost) { 1819e102a522SHong Zhang ghost2 = PETSC_TRUE; 1820e102a522SHong Zhang } else { 1821e102a522SHong Zhang ghost2 = PETSC_FALSE; 1822e102a522SHong Zhang } 1823e102a522SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr); 1824840c2264SHong Zhang } 1825840c2264SHong Zhang 182689898e50SHong Zhang /* Set preallocation for vertex self */ 1827840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1828840c2264SHong Zhang if (!ghost) { 1829840c2264SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 18308675203cSHong Zhang if (network->Jv) { 1831840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 18328675203cSHong Zhang } else Juser = NULL; 1833bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 1834840c2264SHong Zhang } 1835bdcb62a2SHong Zhang if (ghost) { 1836bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1837bdcb62a2SHong Zhang } 1838840c2264SHong Zhang } 1839840c2264SHong Zhang 1840840c2264SHong Zhang ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr); 1841840c2264SHong Zhang ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr); 18425cf7da58SHong Zhang 18435cf7da58SHong Zhang ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr); 18445cf7da58SHong Zhang 18455cf7da58SHong Zhang ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr); 1846840c2264SHong Zhang ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr); 1847840c2264SHong Zhang 1848840c2264SHong Zhang ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr); 1849840c2264SHong Zhang ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr); 1850840c2264SHong Zhang for (j=0; j<localSize; j++) { 1851e102a522SHong Zhang dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 1852e102a522SHong Zhang onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 1853840c2264SHong Zhang } 1854840c2264SHong Zhang ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr); 1855840c2264SHong Zhang ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr); 1856840c2264SHong Zhang ierr = VecDestroy(&vd_nz);CHKERRQ(ierr); 1857840c2264SHong Zhang ierr = VecDestroy(&vo_nz);CHKERRQ(ierr); 1858840c2264SHong Zhang 18595cf7da58SHong Zhang ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr); 18605cf7da58SHong Zhang ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr); 18615cf7da58SHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 18625cf7da58SHong Zhang 18635cf7da58SHong Zhang ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr); 18645cf7da58SHong Zhang 186589898e50SHong Zhang /* (2) Set matrix entries for edges */ 186689898e50SHong Zhang /*----------------------------------*/ 18671ad426b7SHong Zhang for (e=eStart; e<eEnd; e++) { 1868bfbc38dcSHong Zhang /* Get row indices */ 18691ad426b7SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,e,&rstart);CHKERRQ(ierr); 187017df6e9eSHong Zhang ierr = DMNetworkGetNumVariables(dm,e,&nrows);CHKERRQ(ierr); 18714b976069SHong Zhang if (nrows) { 187217df6e9eSHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 18731ad426b7SHong Zhang 1874bfbc38dcSHong Zhang /* Set matrix entries for conntected vertices */ 1875d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 1876bfbc38dcSHong Zhang for (v=0; v<2; v++) { 1877bfbc38dcSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,cone[v],&cstart);CHKERRQ(ierr); 1878883e35e8SHong Zhang ierr = DMNetworkGetNumVariables(dm,cone[v],&ncols);CHKERRQ(ierr); 18793e97b6e8SHong Zhang 18808675203cSHong Zhang if (network->Je) { 1881a4e85ca8SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 18828675203cSHong Zhang } else Juser = NULL; 1883a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 1884bfbc38dcSHong Zhang } 188517df6e9eSHong Zhang 1886bfbc38dcSHong Zhang /* Set matrix entries for edge self */ 18873e97b6e8SHong Zhang cstart = rstart; 18888675203cSHong Zhang if (network->Je) { 1889a4e85ca8SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 18908675203cSHong Zhang } else Juser = NULL; 1891a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 18921ad426b7SHong Zhang } 18934b976069SHong Zhang } 18941ad426b7SHong Zhang 1895bfbc38dcSHong Zhang /* Set matrix entries for vertices */ 189683b2e829SHong Zhang /*---------------------------------*/ 18971ad426b7SHong Zhang for (v=vStart; v<vEnd; v++) { 1898bfbc38dcSHong Zhang /* Get row indices */ 1899596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&rstart);CHKERRQ(ierr); 1900596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,v,&nrows);CHKERRQ(ierr); 19014b976069SHong Zhang if (!nrows) continue; 1902596e729fSHong Zhang 1903bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 1904bdcb62a2SHong Zhang if (ghost) { 1905bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 1906bdcb62a2SHong Zhang } else { 1907bdcb62a2SHong Zhang rows_v = rows; 1908bdcb62a2SHong Zhang } 1909bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 1910596e729fSHong Zhang 1911bfbc38dcSHong Zhang /* Get supporting edges and connected vertices */ 1912596e729fSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1913596e729fSHong Zhang 1914596e729fSHong Zhang for (e=0; e<nedges; e++) { 1915bfbc38dcSHong Zhang /* Supporting edges */ 1916596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,edges[e],&cstart);CHKERRQ(ierr); 1917596e729fSHong Zhang ierr = DMNetworkGetNumVariables(dm,edges[e],&ncols);CHKERRQ(ierr); 1918596e729fSHong Zhang 19198675203cSHong Zhang if (network->Jv) { 1920a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 19218675203cSHong Zhang } else Juser = NULL; 1922bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1923596e729fSHong Zhang 1924bfbc38dcSHong Zhang /* Connected vertices */ 1925d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 19262a945128SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 19272a945128SHong Zhang 192844aca652SHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,vc,&cstart);CHKERRQ(ierr); 192944aca652SHong Zhang ierr = DMNetworkGetNumVariables(dm,vc,&ncols);CHKERRQ(ierr); 1930a4e85ca8SHong Zhang 19318675203cSHong Zhang if (network->Jv) { 1932a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 19338675203cSHong Zhang } else Juser = NULL; 1934bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 1935596e729fSHong Zhang } 1936596e729fSHong Zhang 1937bfbc38dcSHong Zhang /* Set matrix entries for vertex self */ 19381ad426b7SHong Zhang if (!ghost) { 1939596e729fSHong Zhang ierr = DMNetworkGetVariableGlobalOffset(dm,v,&cstart);CHKERRQ(ierr); 19408675203cSHong Zhang if (network->Jv) { 1941a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 19428675203cSHong Zhang } else Juser = NULL; 1943bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr); 1944bdcb62a2SHong Zhang } 1945bdcb62a2SHong Zhang if (ghost) { 1946bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 1947bdcb62a2SHong Zhang } 19481ad426b7SHong Zhang } 1949a4e85ca8SHong Zhang ierr = PetscFree(rows);CHKERRQ(ierr); 1950bdcb62a2SHong Zhang 19511ad426b7SHong Zhang ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 19521ad426b7SHong Zhang ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 1953dd6f46cdSHong Zhang 19545f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 19555f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19565f2c45f1SShri Abhyankar } 19575f2c45f1SShri Abhyankar 19585f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 19595f2c45f1SShri Abhyankar { 19605f2c45f1SShri Abhyankar PetscErrorCode ierr; 19615f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 19622727e31bSShri Abhyankar PetscInt j; 19635f2c45f1SShri Abhyankar 19645f2c45f1SShri Abhyankar PetscFunctionBegin; 19658415c774SShri Abhyankar if (--network->refct > 0) PetscFunctionReturn(0); 196683b2e829SHong Zhang if (network->Je) { 196783b2e829SHong Zhang ierr = PetscFree(network->Je);CHKERRQ(ierr); 196883b2e829SHong Zhang } 196983b2e829SHong Zhang if (network->Jv) { 1970883e35e8SHong Zhang ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 197183b2e829SHong Zhang ierr = PetscFree(network->Jv);CHKERRQ(ierr); 19721ad426b7SHong Zhang } 197313c2a604SAdrian Maldonado 197413c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr); 197513c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr); 197613c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr); 197713c2a604SAdrian Maldonado if (network->vertex.sf) { 197813c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr); 197913c2a604SAdrian Maldonado } 198013c2a604SAdrian Maldonado /* edge */ 198113c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr); 198213c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr); 198313c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr); 198413c2a604SAdrian Maldonado if (network->edge.sf) { 198513c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr); 198613c2a604SAdrian Maldonado } 19875f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 19885f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 19895f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 199083b2e829SHong Zhang 19912727e31bSShri Abhyankar for(j=0; j<network->nsubnet; j++) { 19922727e31bSShri Abhyankar ierr = PetscFree(network->subnet[j].edges);CHKERRQ(ierr); 19932727e31bSShri Abhyankar } 1994caf410d2SHong Zhang ierr = PetscFree(network->subnetvtx);CHKERRQ(ierr); 1995caf410d2SHong Zhang 1996e2aaf10cSShri Abhyankar ierr = PetscFree(network->subnet);CHKERRQ(ierr); 19975f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 1998caf410d2SHong Zhang ierr = PetscFree2(network->header,network->cvalue);CHKERRQ(ierr); 19995f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 20005f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20015f2c45f1SShri Abhyankar } 20025f2c45f1SShri Abhyankar 20035f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm,PetscViewer viewer) 20045f2c45f1SShri Abhyankar { 20055f2c45f1SShri Abhyankar PetscErrorCode ierr; 20065f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 2007caf410d2SHong Zhang PetscBool iascii; 2008caf410d2SHong Zhang PetscMPIInt rank; 2009caf410d2SHong Zhang PetscInt p,nsubnet; 20105f2c45f1SShri Abhyankar 20115f2c45f1SShri Abhyankar PetscFunctionBegin; 2012caf410d2SHong Zhang if (!dm->setupcalled) SETERRQ(PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first"); 2013caf410d2SHong Zhang ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRQ(ierr); 2014caf410d2SHong Zhang PetscValidHeaderSpecific(dm,DM_CLASSID, 1); 2015caf410d2SHong Zhang PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2016caf410d2SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 2017caf410d2SHong Zhang if (iascii) { 2018caf410d2SHong Zhang const PetscInt *cone,*vtx,*edges; 2019caf410d2SHong Zhang PetscInt vfrom,vto,i,j,nv,ne; 2020caf410d2SHong Zhang 2021caf410d2SHong Zhang nsubnet = network->nsubnet - network->ncsubnet; /* num of subnetworks */ 2022caf410d2SHong Zhang ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 2023caf410d2SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " [%d] nsubnet: %D; nsubnetCouple: %D; nEdges: %D; nVertices: %D\n",rank,nsubnet,network->ncsubnet,network->nEdges,network->nVertices);CHKERRQ(ierr); 2024caf410d2SHong Zhang 2025caf410d2SHong Zhang for (i=0; i<nsubnet; i++) { 2026caf410d2SHong Zhang ierr = DMNetworkGetSubnetworkInfo(dm,i,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 2027caf410d2SHong Zhang if (ne) { 2028caf410d2SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %D: nEdges %D, nVertices %D\n",i,ne,nv);CHKERRQ(ierr); 2029caf410d2SHong Zhang for (j=0; j<ne; j++) { 2030caf410d2SHong Zhang p = edges[j]; 2031caf410d2SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,p,&cone);CHKERRQ(ierr); 2032caf410d2SHong Zhang ierr = DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);CHKERRQ(ierr); 2033caf410d2SHong Zhang ierr = DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);CHKERRQ(ierr); 2034caf410d2SHong Zhang ierr = DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p);CHKERRQ(ierr); 2035caf410d2SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " edge %D: %D----> %D\n",p,vfrom,vto);CHKERRQ(ierr); 2036caf410d2SHong Zhang } 2037caf410d2SHong Zhang } 2038caf410d2SHong Zhang } 2039caf410d2SHong Zhang /* Coupling subnets */ 2040caf410d2SHong Zhang nsubnet = network->nsubnet; 2041caf410d2SHong Zhang for (; i<nsubnet; i++) { 2042caf410d2SHong Zhang ierr = DMNetworkGetSubnetworkInfo(dm,i,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 2043caf410d2SHong Zhang if (ne) { 2044caf410d2SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %D (couple): nEdges %D, nVertices %D\n",i,ne,nv);CHKERRQ(ierr); 2045caf410d2SHong Zhang for (j=0; j<ne; j++) { 2046caf410d2SHong Zhang p = edges[j]; 2047caf410d2SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,p,&cone);CHKERRQ(ierr); 2048caf410d2SHong Zhang ierr = DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);CHKERRQ(ierr); 2049caf410d2SHong Zhang ierr = DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);CHKERRQ(ierr); 2050caf410d2SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " edge %D: %D----> %D\n",p,vfrom,vto);CHKERRQ(ierr); 2051caf410d2SHong Zhang } 2052caf410d2SHong Zhang } 2053caf410d2SHong Zhang } 2054caf410d2SHong Zhang ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2055caf410d2SHong Zhang ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 2056caf410d2SHong Zhang } else SETERRQ1(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name); 20575f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20585f2c45f1SShri Abhyankar } 20595f2c45f1SShri Abhyankar 20605f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 20615f2c45f1SShri Abhyankar { 20625f2c45f1SShri Abhyankar PetscErrorCode ierr; 20635f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 20645f2c45f1SShri Abhyankar 20655f2c45f1SShri Abhyankar PetscFunctionBegin; 20665f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 20675f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20685f2c45f1SShri Abhyankar } 20695f2c45f1SShri Abhyankar 20705f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 20715f2c45f1SShri Abhyankar { 20725f2c45f1SShri Abhyankar PetscErrorCode ierr; 20735f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 20745f2c45f1SShri Abhyankar 20755f2c45f1SShri Abhyankar PetscFunctionBegin; 20765f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 20775f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20785f2c45f1SShri Abhyankar } 20795f2c45f1SShri Abhyankar 20805f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 20815f2c45f1SShri Abhyankar { 20825f2c45f1SShri Abhyankar PetscErrorCode ierr; 20835f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 20845f2c45f1SShri Abhyankar 20855f2c45f1SShri Abhyankar PetscFunctionBegin; 20865f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 20875f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20885f2c45f1SShri Abhyankar } 20895f2c45f1SShri Abhyankar 20905f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 20915f2c45f1SShri Abhyankar { 20925f2c45f1SShri Abhyankar PetscErrorCode ierr; 20935f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 20945f2c45f1SShri Abhyankar 20955f2c45f1SShri Abhyankar PetscFunctionBegin; 20965f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 20975f2c45f1SShri Abhyankar PetscFunctionReturn(0); 20985f2c45f1SShri Abhyankar } 2099