1af0996ceSBarry Smith #include <petsc/private/dmnetworkimpl.h> /*I "petscdmnetwork.h" I*/ 25f2c45f1SShri Abhyankar 354dfd506SHong Zhang /* 454dfd506SHong Zhang Creates the component header and value objects for a network point 554dfd506SHong Zhang */ 654dfd506SHong Zhang static PetscErrorCode SetUpNetworkHeaderComponentValue(DM dm,DMNetworkComponentHeader header,DMNetworkComponentValue cvalue) 754dfd506SHong Zhang { 854dfd506SHong Zhang PetscErrorCode ierr; 954dfd506SHong Zhang 1054dfd506SHong Zhang PetscFunctionBegin; 1154dfd506SHong Zhang /* Allocate arrays for component information */ 1254dfd506SHong Zhang ierr = PetscCalloc5(header->maxcomps,&header->size,header->maxcomps,&header->key,header->maxcomps,&header->offset,header->maxcomps,&header->nvar,header->maxcomps,&header->offsetvarrel);CHKERRQ(ierr); 1354dfd506SHong Zhang ierr = PetscCalloc1(header->maxcomps,&cvalue->data);CHKERRQ(ierr); 1454dfd506SHong Zhang 1554dfd506SHong Zhang /* The size of the header is the size of struct _p_DMNetworkComponentHeader. Since the struct contains PetscInt pointers we cannot use sizeof(struct). So, we need to explicitly calculate the size. 1654dfd506SHong Zhang If the data header struct changes then this header size calculation needs to be updated. */ 1754dfd506SHong Zhang header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt); 1854dfd506SHong Zhang header->hsize /= sizeof(DMNetworkComponentGenericDataType); 1954dfd506SHong Zhang PetscFunctionReturn(0); 2054dfd506SHong Zhang } 2154dfd506SHong Zhang 225f2c45f1SShri Abhyankar /*@ 23556ed216SShri Abhyankar DMNetworkGetPlex - Gets the Plex DM associated with this network DM 24556ed216SShri Abhyankar 25556ed216SShri Abhyankar Not collective 26556ed216SShri Abhyankar 27556ed216SShri Abhyankar Input Parameters: 282bf73ac6SHong Zhang . dm - the dm object 292bf73ac6SHong Zhang 302bf73ac6SHong Zhang Output Parameters: 312bf73ac6SHong Zhang . plexdm - the plex dm object 32556ed216SShri Abhyankar 33556ed216SShri Abhyankar Level: Advanced 34556ed216SShri Abhyankar 35556ed216SShri Abhyankar .seealso: DMNetworkCreate() 36556ed216SShri Abhyankar @*/ 372bf73ac6SHong Zhang PetscErrorCode DMNetworkGetPlex(DM dm,DM *plexdm) 38556ed216SShri Abhyankar { 392bf73ac6SHong Zhang DM_Network *network = (DM_Network*)dm->data; 40556ed216SShri Abhyankar 41556ed216SShri Abhyankar PetscFunctionBegin; 42556ed216SShri Abhyankar *plexdm = network->plex; 43556ed216SShri Abhyankar PetscFunctionReturn(0); 44556ed216SShri Abhyankar } 45556ed216SShri Abhyankar 46556ed216SShri Abhyankar /*@ 472bf73ac6SHong Zhang DMNetworkGetNumSubNetworks - Gets the the number of subnetworks 4872c3e9f7SHong Zhang 492bf73ac6SHong Zhang Not collective 5072c3e9f7SHong Zhang 51f899ff85SJose E. Roman Input Parameter: 522bf73ac6SHong Zhang . dm - the dm object 532bf73ac6SHong Zhang 542bf73ac6SHong Zhang Output Parameters: 552bf73ac6SHong Zhang + nsubnet - local number of subnetworks 562bf73ac6SHong Zhang - Nsubnet - global number of subnetworks 5772c3e9f7SHong Zhang 5897bb938eSShri Abhyankar Level: beginner 5972c3e9f7SHong Zhang 602bf73ac6SHong Zhang .seealso: DMNetworkCreate(), DMNetworkSetNumSubNetworks() 6172c3e9f7SHong Zhang @*/ 622bf73ac6SHong Zhang PetscErrorCode DMNetworkGetNumSubNetworks(DM dm,PetscInt *nsubnet,PetscInt *Nsubnet) 6372c3e9f7SHong Zhang { 642bf73ac6SHong Zhang DM_Network *network = (DM_Network*)dm->data; 6572c3e9f7SHong Zhang 6672c3e9f7SHong Zhang PetscFunctionBegin; 672bf73ac6SHong Zhang if (nsubnet) *nsubnet = network->nsubnet; 682bf73ac6SHong Zhang if (Nsubnet) *Nsubnet = network->Nsubnet; 6972c3e9f7SHong Zhang PetscFunctionReturn(0); 7072c3e9f7SHong Zhang } 7172c3e9f7SHong Zhang 7272c3e9f7SHong Zhang /*@ 732bf73ac6SHong Zhang DMNetworkSetNumSubNetworks - Sets the number of subnetworks 745f2c45f1SShri Abhyankar 75d083f849SBarry Smith Collective on dm 765f2c45f1SShri Abhyankar 775f2c45f1SShri Abhyankar Input Parameters: 785f2c45f1SShri Abhyankar + dm - the dm object 792bf73ac6SHong Zhang . nsubnet - local number of subnetworks 802bf73ac6SHong Zhang - Nsubnet - global number of subnetworks 815f2c45f1SShri Abhyankar 8297bb938eSShri Abhyankar Level: beginner 831b266c99SBarry Smith 842bf73ac6SHong Zhang .seealso: DMNetworkCreate(), DMNetworkGetNumSubNetworks() 855f2c45f1SShri Abhyankar @*/ 862bf73ac6SHong Zhang PetscErrorCode DMNetworkSetNumSubNetworks(DM dm,PetscInt nsubnet,PetscInt Nsubnet) 875f2c45f1SShri Abhyankar { 885f2c45f1SShri Abhyankar PetscErrorCode ierr; 895f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 905f2c45f1SShri Abhyankar 915f2c45f1SShri Abhyankar PetscFunctionBegin; 929ace16cdSJacob Faibussowitsch PetscAssertFalse(network->Nsubnet != 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Network sizes alread set, cannot resize the network"); 932bf73ac6SHong Zhang 945f2c45f1SShri Abhyankar PetscValidHeaderSpecific(dm,DM_CLASSID,1); 952bf73ac6SHong Zhang PetscValidLogicalCollectiveInt(dm,nsubnet,2); 962bf73ac6SHong Zhang PetscValidLogicalCollectiveInt(dm,Nsubnet,3); 977765340cSHong Zhang 982bf73ac6SHong Zhang if (Nsubnet == PETSC_DECIDE) { 999ace16cdSJacob Faibussowitsch PetscAssertFalse(nsubnet < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Number of local subnetworks %D cannot be less than 0",nsubnet); 100820f2d46SBarry Smith ierr = MPIU_Allreduce(&nsubnet,&Nsubnet,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr); 101caf410d2SHong Zhang } 1029ace16cdSJacob Faibussowitsch PetscAssertFalse(Nsubnet < 1,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_INCOMP,"Number of global subnetworks %D cannot be less than 1",Nsubnet); 103caf410d2SHong Zhang 1042bf73ac6SHong Zhang network->Nsubnet = Nsubnet; 1052bf73ac6SHong Zhang network->nsubnet = 0; /* initia value; will be determind by DMNetworkAddSubnetwork() */ 1062bf73ac6SHong Zhang ierr = PetscCalloc1(Nsubnet,&network->subnet);CHKERRQ(ierr); 107caf410d2SHong Zhang 1082bf73ac6SHong Zhang /* num of shared vertices */ 1092bf73ac6SHong Zhang network->nsvtx = 0; 1102bf73ac6SHong Zhang network->Nsvtx = 0; 1117765340cSHong Zhang PetscFunctionReturn(0); 1127765340cSHong Zhang } 1137765340cSHong Zhang 1145f2c45f1SShri Abhyankar /*@ 1152bf73ac6SHong Zhang DMNetworkAddSubnetwork - Add a subnetwork 1165f2c45f1SShri Abhyankar 1172bf73ac6SHong Zhang Collective on dm 1185f2c45f1SShri Abhyankar 1195f2c45f1SShri Abhyankar Input Parameters: 120e3e68989SHong Zhang + dm - the dm object 1212bf73ac6SHong Zhang . name - name of the subnetwork 1222bf73ac6SHong Zhang . ne - number of local edges of this subnetwork 1238d1cd658SBarry Smith - edgelist - list of edges for this subnetwork, this is a one dimensional array with pairs of entries being the two vertices (in global numbering of the vertices) of each edge, 1248d1cd658SBarry Smith $ [first vertex of first edge, second vertex of first edge, first vertex of second edge, second vertex of second edge, etc] 1252bf73ac6SHong Zhang 1262bf73ac6SHong Zhang Output Parameters: 1272bf73ac6SHong Zhang . netnum - global index of the subnetwork 1285f2c45f1SShri Abhyankar 1295f2c45f1SShri Abhyankar Notes: 1305f2c45f1SShri Abhyankar There is no copy involved in this operation, only the pointer is referenced. The edgelist should 1312bf73ac6SHong Zhang not be destroyed before the call to DMNetworkLayoutSetUp() 1325f2c45f1SShri Abhyankar 1338d1cd658SBarry Smith A network can comprise of a single subnetwork OR multiple subnetworks. For a single subnetwork, the subnetwork can be read either in serial or parallel. For a multiple subnetworks, 1348d1cd658SBarry Smith each subnetwork topology needs to be set on a unique rank and the communicator size needs to be at least equal to the number of subnetworks. 135f11a936eSBarry Smith 13697bb938eSShri Abhyankar Level: beginner 1375f2c45f1SShri Abhyankar 138e3e68989SHong Zhang Example usage: 139f11a936eSBarry Smith Consider the following networks: 140f11a936eSBarry Smith 1) A sigle subnetwork: 141e3e68989SHong Zhang .vb 142f11a936eSBarry Smith network 0: 143f11a936eSBarry Smith rank[0]: 144f11a936eSBarry Smith v0 --> v2; v1 --> v2 145f11a936eSBarry Smith rank[1]: 146f11a936eSBarry Smith v3 --> v5; v4 --> v5 147e3e68989SHong Zhang .ve 148e3e68989SHong Zhang 149e3e68989SHong Zhang The resulting input 150f11a936eSBarry Smith network 0: 151f11a936eSBarry Smith rank[0]: 152f11a936eSBarry Smith ne = 2 153f11a936eSBarry Smith edgelist = [0 2 | 1 2] 154f11a936eSBarry Smith rank[1]: 155f11a936eSBarry Smith ne = 2 156f11a936eSBarry Smith edgelist = [3 5 | 4 5] 157f11a936eSBarry Smith 158f11a936eSBarry Smith 2) Two subnetworks: 159f11a936eSBarry Smith .vb 160f11a936eSBarry Smith subnetwork 0: 161f11a936eSBarry Smith rank[0]: 162f11a936eSBarry Smith v0 --> v2; v2 --> v1; v1 --> v3; 163f11a936eSBarry Smith subnetwork 1: 164f11a936eSBarry Smith rank[1]: 165f11a936eSBarry Smith v0 --> v3; v3 --> v2; v2 --> v1; 166f11a936eSBarry Smith .ve 167f11a936eSBarry Smith 168f11a936eSBarry Smith The resulting input 169f11a936eSBarry Smith subnetwork 0: 170f11a936eSBarry Smith rank[0]: 171f11a936eSBarry Smith ne = 3 172f11a936eSBarry Smith edgelist = [0 2 | 2 1 | 1 3] 173f11a936eSBarry Smith rank[1]: 174f11a936eSBarry Smith ne = 0 175f11a936eSBarry Smith edgelist = NULL 176f11a936eSBarry Smith 177f11a936eSBarry Smith subnetwork 1: 178f11a936eSBarry Smith rank[0]: 179f11a936eSBarry Smith ne = 0 180f11a936eSBarry Smith edgelist = NULL 181f11a936eSBarry Smith rank[1]: 182f11a936eSBarry Smith edgelist = [0 3 | 3 2 | 2 1] 183e3e68989SHong Zhang 1842bf73ac6SHong Zhang .seealso: DMNetworkCreate(), DMNetworkSetNumSubnetworks() 1855f2c45f1SShri Abhyankar @*/ 186f11a936eSBarry Smith PetscErrorCode DMNetworkAddSubnetwork(DM dm,const char* name,PetscInt ne,PetscInt edgelist[],PetscInt *netnum) 1875f2c45f1SShri Abhyankar { 1882bf73ac6SHong Zhang PetscErrorCode ierr; 1895f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 190f11a936eSBarry Smith PetscInt i,Nedge,j,Nvtx,nvtx; 191f11a936eSBarry Smith PetscBT table; 1925f2c45f1SShri Abhyankar 1935f2c45f1SShri Abhyankar PetscFunctionBegin; 1948d1cd658SBarry Smith for (i=0; i<ne; i++) { 1959ace16cdSJacob Faibussowitsch PetscAssertFalse(edgelist[2*i] == edgelist[2*i+1],PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"Edge %D has the same vertex %D at each endpoint",i,edgelist[2*i]); 1968d1cd658SBarry Smith } 197f11a936eSBarry Smith /* Get global total Nvtx = max(edgelist[])+1 for this subnet */ 198f11a936eSBarry Smith nvtx = -1; i = 0; 199f11a936eSBarry Smith for (j=0; j<ne; j++) { 200f11a936eSBarry Smith nvtx = PetscMax(nvtx, edgelist[i]); i++; 201f11a936eSBarry Smith nvtx = PetscMax(nvtx, edgelist[i]); i++; 202f11a936eSBarry Smith } 203f11a936eSBarry Smith nvtx++; 204f11a936eSBarry Smith ierr = MPIU_Allreduce(&nvtx,&Nvtx,1,MPIU_INT,MPI_MAX,PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr); 205f11a936eSBarry Smith 206f11a936eSBarry Smith /* Get local nvtx for this subnet */ 207f11a936eSBarry Smith ierr = PetscBTCreate(Nvtx,&table);CHKERRQ(ierr); 208f11a936eSBarry Smith ierr = PetscBTMemzero(Nvtx,table);CHKERRQ(ierr); 209f11a936eSBarry Smith i = 0; 210f11a936eSBarry Smith for (j=0; j<ne; j++) { 211f11a936eSBarry Smith ierr = PetscBTSet(table,edgelist[i]);CHKERRQ(ierr); 212f11a936eSBarry Smith i++; 213f11a936eSBarry Smith ierr = PetscBTSet(table,edgelist[i]);CHKERRQ(ierr); 214f11a936eSBarry Smith i++; 215f11a936eSBarry Smith } 216f11a936eSBarry Smith nvtx = 0; 217f11a936eSBarry Smith for (j=0; j<Nvtx; j++) { 218f11a936eSBarry Smith if (PetscBTLookup(table,j)) nvtx++; 219f11a936eSBarry Smith } 220f11a936eSBarry Smith ierr = PetscBTDestroy(&table);CHKERRQ(ierr); 221f11a936eSBarry Smith 222f11a936eSBarry Smith /* Get global total Nedge for this subnet */ 223f11a936eSBarry Smith ierr = MPIU_Allreduce(&ne,&Nedge,1,MPIU_INT,MPI_SUM,PetscObjectComm((PetscObject)dm));CHKERRMPI(ierr); 224f11a936eSBarry Smith 225f11a936eSBarry Smith i = network->nsubnet; 2262bf73ac6SHong Zhang if (name) { 2272bf73ac6SHong Zhang ierr = PetscStrcpy(network->subnet[i].name,name);CHKERRQ(ierr); 228e3e68989SHong Zhang } 229f11a936eSBarry Smith network->subnet[i].nvtx = nvtx; /* include ghost vertices */ 2302bf73ac6SHong Zhang network->subnet[i].nedge = ne; 2312bf73ac6SHong Zhang network->subnet[i].edgelist = edgelist; 232f11a936eSBarry Smith network->subnet[i].Nvtx = Nvtx; 233f11a936eSBarry Smith network->subnet[i].Nedge = Nedge; 2342bf73ac6SHong Zhang 2352bf73ac6SHong Zhang /* ---------------------------------------------------------- 2362bf73ac6SHong Zhang p=v or e; 2372bf73ac6SHong Zhang subnet[0].pStart = 0 2382bf73ac6SHong Zhang subnet[i+1].pStart = subnet[i].pEnd = subnet[i].pStart + (nE[i] or NV[i]) 2392bf73ac6SHong Zhang ----------------------------------------------------------------------- */ 2402bf73ac6SHong Zhang /* GLOBAL subnet[].vStart and vEnd, used by DMNetworkLayoutSetUp() */ 2412bf73ac6SHong Zhang network->subnet[i].vStart = network->NVertices; 2422bf73ac6SHong Zhang network->subnet[i].vEnd = network->subnet[i].vStart + network->subnet[i].Nvtx; /* global vEnd of subnet[i] */ 2432bf73ac6SHong Zhang 244f11a936eSBarry Smith network->nVertices += nvtx; /* include ghost vertices */ 2452bf73ac6SHong Zhang network->NVertices += network->subnet[i].Nvtx; 2462bf73ac6SHong Zhang 2472bf73ac6SHong Zhang /* LOCAL subnet[].eStart and eEnd, used by DMNetworkLayoutSetUp() */ 2482bf73ac6SHong Zhang network->subnet[i].eStart = network->nEdges; 2492bf73ac6SHong Zhang network->subnet[i].eEnd = network->subnet[i].eStart + ne; 2502bf73ac6SHong Zhang network->nEdges += ne; 2512bf73ac6SHong Zhang network->NEdges += network->subnet[i].Nedge; 2522bf73ac6SHong Zhang 2532bf73ac6SHong Zhang ierr = PetscStrcpy(network->subnet[i].name,name);CHKERRQ(ierr); 2542bf73ac6SHong Zhang if (netnum) *netnum = network->nsubnet; 2552bf73ac6SHong Zhang network->nsubnet++; 2562bf73ac6SHong Zhang PetscFunctionReturn(0); 2572bf73ac6SHong Zhang } 2582bf73ac6SHong Zhang 2592bf73ac6SHong Zhang /* 2602bf73ac6SHong Zhang SetUp a single svtx struct. See SVtx defined in dmnetworkimpl.h 2612bf73ac6SHong Zhang Set gidx and type if input v=(net,idx) is a from_vertex; 2622bf73ac6SHong Zhang Get gid, type and index in the svtx array if input v=(net,idx) is a to_vertex. 2632bf73ac6SHong Zhang 2642bf73ac6SHong Zhang Input: Nsvtx, svtx, net, idx, gidx 2652bf73ac6SHong Zhang Output: gidx, svtype, svtx_idx 2662bf73ac6SHong Zhang */ 2672bf73ac6SHong Zhang static PetscErrorCode SVtxSetUp(PetscInt Nsvtx,SVtx *svtx,PetscInt net,PetscInt idx,PetscInt *gidx,SVtxType *svtype,PetscInt *svtx_idx) 2682bf73ac6SHong Zhang { 2692bf73ac6SHong Zhang PetscInt i,j,*svto; 2702bf73ac6SHong Zhang SVtxType vtype; 2712bf73ac6SHong Zhang 2722bf73ac6SHong Zhang PetscFunctionBegin; 2732bf73ac6SHong Zhang if (!Nsvtx) PetscFunctionReturn(0); 2742bf73ac6SHong Zhang 2752bf73ac6SHong Zhang vtype = SVNONE; 2762bf73ac6SHong Zhang for (i=0; i<Nsvtx; i++) { 2772bf73ac6SHong Zhang if (net == svtx[i].sv[0] && idx == svtx[i].sv[1]) { 2782bf73ac6SHong Zhang /* (1) input vertex net.idx is a shared from_vertex, set its global index and output its svtype */ 2792bf73ac6SHong Zhang svtx[i].gidx = *gidx; /* set gidx */ 2802bf73ac6SHong Zhang vtype = SVFROM; 2812bf73ac6SHong Zhang } else { /* loop over svtx[i].n */ 2822bf73ac6SHong Zhang for (j=1; j<svtx[i].n; j++) { 2832bf73ac6SHong Zhang svto = svtx[i].sv + 2*j; 2842bf73ac6SHong Zhang if (net == svto[0] && idx == svto[1]) { 2852bf73ac6SHong Zhang /* input vertex net.idx is a shared to_vertex, output its global index and its svtype */ 2862bf73ac6SHong Zhang *gidx = svtx[i].gidx; /* output gidx for to_vertex */ 2872bf73ac6SHong Zhang vtype = SVTO; 2882bf73ac6SHong Zhang } 2892bf73ac6SHong Zhang } 2902bf73ac6SHong Zhang } 2912bf73ac6SHong Zhang if (vtype != SVNONE) break; 2922bf73ac6SHong Zhang } 2932bf73ac6SHong Zhang if (svtype) *svtype = vtype; 2942bf73ac6SHong Zhang if (svtx_idx) *svtx_idx = i; 2952bf73ac6SHong Zhang PetscFunctionReturn(0); 2962bf73ac6SHong Zhang } 2972bf73ac6SHong Zhang 2982bf73ac6SHong Zhang /* 2992bf73ac6SHong Zhang Add a new shared vertex sv=(net,idx) to table svtas[ita] 3002bf73ac6SHong Zhang */ 3012bf73ac6SHong Zhang static PetscErrorCode TableAddSVtx(PetscTable *svtas,PetscInt ita,PetscInt* tdata,PetscInt *sv_wk,PetscInt *ii,PetscInt *sedgelist,PetscInt k,DM_Network *network,PetscInt **ta2sv) 3022bf73ac6SHong Zhang { 3032bf73ac6SHong Zhang PetscInt net,idx,gidx,i=*ii; 3042bf73ac6SHong Zhang PetscErrorCode ierr; 3052bf73ac6SHong Zhang 3062bf73ac6SHong Zhang PetscFunctionBegin; 3072bf73ac6SHong Zhang net = sv_wk[2*i] = sedgelist[k]; 3082bf73ac6SHong Zhang idx = sv_wk[2*i+1] = sedgelist[k+1]; 3092bf73ac6SHong Zhang gidx = network->subnet[net].vStart + idx; 3102bf73ac6SHong Zhang ierr = PetscTableAdd(svtas[ita],gidx+1,tdata[ita]+1,INSERT_VALUES);CHKERRQ(ierr); 3112bf73ac6SHong Zhang *(ta2sv[ita] + tdata[ita]) = i; /* maps tdata to index of sv_wk; sv_wk keeps (net,idx) info */ 3122bf73ac6SHong Zhang tdata[ita]++; (*ii)++; 3132bf73ac6SHong Zhang PetscFunctionReturn(0); 3142bf73ac6SHong Zhang } 3152bf73ac6SHong Zhang 3162bf73ac6SHong Zhang /* 3172bf73ac6SHong Zhang Create an array of shared vertices. See SVtx and SVtxType in dmnetworkimpl.h 3182bf73ac6SHong Zhang 3192bf73ac6SHong Zhang Input: dm, Nsedgelist, sedgelist 3202bf73ac6SHong Zhang Output: Nsvtx,svtx 3212bf73ac6SHong Zhang 3222bf73ac6SHong Zhang Note: Output svtx is organized as 3232bf73ac6SHong Zhang sv(net[0],idx[0]) --> sv(net[1],idx[1]) 3242bf73ac6SHong Zhang --> sv(net[1],idx[1]) 3252bf73ac6SHong Zhang ... 3262bf73ac6SHong Zhang --> sv(net[n-1],idx[n-1]) 3272bf73ac6SHong Zhang and net[0] < net[1] < ... < net[n-1] 3282bf73ac6SHong Zhang where sv[0] has SVFROM type, sv[i], i>0, has SVTO type. 3292bf73ac6SHong Zhang */ 3302bf73ac6SHong Zhang static PetscErrorCode SVtxCreate(DM dm,PetscInt Nsedgelist,PetscInt *sedgelist,PetscInt *Nsvtx,SVtx **svtx) 3312bf73ac6SHong Zhang { 3322bf73ac6SHong Zhang PetscErrorCode ierr; 3332bf73ac6SHong Zhang SVtx *sedges = NULL; 3342bf73ac6SHong Zhang PetscInt *sv,k,j,nsv,*tdata,**ta2sv; 3352bf73ac6SHong Zhang PetscTable *svtas; 3362bf73ac6SHong Zhang PetscInt gidx,net,idx,i,nta,ita,idx_from,idx_to,n,*sv_wk; 3372bf73ac6SHong Zhang DM_Network *network = (DM_Network*)dm->data; 3382bf73ac6SHong Zhang PetscTablePosition ppos; 3392bf73ac6SHong Zhang 3402bf73ac6SHong Zhang PetscFunctionBegin; 3412bf73ac6SHong Zhang /* (1) Crete ctables svtas */ 3422bf73ac6SHong Zhang ierr = PetscCalloc4(Nsedgelist,&svtas,Nsedgelist,&tdata,4*Nsedgelist,&sv_wk,2*Nsedgelist,&ta2sv);CHKERRQ(ierr); 3432bf73ac6SHong Zhang 3442bf73ac6SHong Zhang j = 0; /* sedgelist counter */ 3452bf73ac6SHong Zhang k = 0; /* sedgelist vertex counter j = 4*k */ 3462bf73ac6SHong Zhang i = 0; /* sv_wk (vertices added to the ctables) counter */ 3472bf73ac6SHong Zhang nta = 0; /* num of sv tables created */ 3482bf73ac6SHong Zhang 3492bf73ac6SHong Zhang /* for j=0 */ 3502bf73ac6SHong Zhang ierr = PetscTableCreate(2*Nsedgelist,network->NVertices+1,&svtas[nta]);CHKERRQ(ierr); 3512bf73ac6SHong Zhang ierr = PetscMalloc1(2*Nsedgelist,&ta2sv[nta]);CHKERRQ(ierr); 3522bf73ac6SHong Zhang 3532bf73ac6SHong Zhang ierr = TableAddSVtx(svtas,nta,tdata,sv_wk,&i,sedgelist,k,network,ta2sv);CHKERRQ(ierr); 3542bf73ac6SHong Zhang ierr = TableAddSVtx(svtas,nta,tdata,sv_wk,&i,sedgelist,k+2,network,ta2sv);CHKERRQ(ierr); 3552bf73ac6SHong Zhang nta++; k += 4; 3562bf73ac6SHong Zhang 3572bf73ac6SHong Zhang for (j = 1; j < Nsedgelist; j++) { 3582bf73ac6SHong Zhang for (ita = 0; ita < nta; ita++) { 3592bf73ac6SHong Zhang /* vfrom */ 3602bf73ac6SHong Zhang net = sedgelist[k]; idx = sedgelist[k+1]; 3612bf73ac6SHong Zhang gidx = network->subnet[net].vStart + idx; /* global index of the vertex net.idx before merging shared vertices */ 3622bf73ac6SHong Zhang ierr = PetscTableFind(svtas[ita],gidx+1,&idx_from);CHKERRQ(ierr); 3632bf73ac6SHong Zhang 3642bf73ac6SHong Zhang /* vto */ 3652bf73ac6SHong Zhang net = sedgelist[k+2]; idx = sedgelist[k+3]; 3662bf73ac6SHong Zhang gidx = network->subnet[net].vStart + idx; 3672bf73ac6SHong Zhang ierr = PetscTableFind(svtas[ita],gidx+1,&idx_to);CHKERRQ(ierr); 3682bf73ac6SHong Zhang 3692bf73ac6SHong Zhang if (idx_from || idx_to) { /* vfrom or vto is on table svtas[ita] */ 3702bf73ac6SHong Zhang idx_from--; idx_to--; 3712bf73ac6SHong Zhang if (idx_from < 0) { /* vto is on svtas[ita] */ 3722bf73ac6SHong Zhang ierr = TableAddSVtx(svtas,ita,tdata,sv_wk,&i,sedgelist,k,network,ta2sv);CHKERRQ(ierr); 3732bf73ac6SHong Zhang break; 3742bf73ac6SHong Zhang } else if (idx_to < 0) { 3752bf73ac6SHong Zhang ierr = TableAddSVtx(svtas,ita,tdata,sv_wk,&i,sedgelist,k+2,network,ta2sv);CHKERRQ(ierr); 3762bf73ac6SHong Zhang break; 3772bf73ac6SHong Zhang } 3782bf73ac6SHong Zhang } 3792bf73ac6SHong Zhang } 3802bf73ac6SHong Zhang 3812bf73ac6SHong Zhang if (ita == nta) { 3822bf73ac6SHong Zhang ierr = PetscTableCreate(2*Nsedgelist,network->NVertices+1,&svtas[nta]);CHKERRQ(ierr); 3832bf73ac6SHong Zhang ierr = PetscMalloc1(2*Nsedgelist, &ta2sv[nta]);CHKERRQ(ierr); 3842bf73ac6SHong Zhang 3852bf73ac6SHong Zhang ierr = TableAddSVtx(svtas,nta,tdata,sv_wk,&i,sedgelist,k,network,ta2sv);CHKERRQ(ierr); 3862bf73ac6SHong Zhang ierr = TableAddSVtx(svtas,nta,tdata,sv_wk,&i,sedgelist,k+2,network,ta2sv);CHKERRQ(ierr); 3872bf73ac6SHong Zhang nta++; 3882bf73ac6SHong Zhang } 3892bf73ac6SHong Zhang k += 4; 3902bf73ac6SHong Zhang } 3912bf73ac6SHong Zhang 3922bf73ac6SHong Zhang /* (2) Construct sedges from ctable 3932bf73ac6SHong Zhang sedges: edges connect vertex sv[0]=(net[0],idx[0]) to vertices sv[k], k=1,...,n-1; 3944f5c2772SJose E. Roman net[k], k=0, ...,n-1, are in ascending order */ 3952bf73ac6SHong Zhang ierr = PetscMalloc1(nta,&sedges);CHKERRQ(ierr); 3962bf73ac6SHong Zhang for (nsv = 0; nsv < nta; nsv++) { 3974f5c2772SJose E. Roman /* for a single svtx, put shared vertices in ascending order of gidx */ 3984f5c2772SJose E. Roman ierr = PetscTableGetCount(svtas[nsv],&n);CHKERRQ(ierr); 3992bf73ac6SHong Zhang ierr = PetscCalloc1(2*n,&sv);CHKERRQ(ierr); 4002bf73ac6SHong Zhang sedges[nsv].sv = sv; 4012bf73ac6SHong Zhang sedges[nsv].n = n; 4022bf73ac6SHong Zhang sedges[nsv].gidx = -1; /* initialization */ 4032bf73ac6SHong Zhang 4042bf73ac6SHong Zhang ierr = PetscTableGetHeadPosition(svtas[nsv],&ppos);CHKERRQ(ierr); 4054f5c2772SJose E. Roman for (k=0; k<n; k++) { /* gidx is sorted in ascending order */ 4062bf73ac6SHong Zhang ierr = PetscTableGetNext(svtas[nsv],&ppos,&gidx,&i);CHKERRQ(ierr); 4072bf73ac6SHong Zhang gidx--; i--; 4082bf73ac6SHong Zhang 4092bf73ac6SHong Zhang j = ta2sv[nsv][i]; /* maps i to index of sv_wk */ 4102bf73ac6SHong Zhang sv[2*k] = sv_wk[2*j]; 4112bf73ac6SHong Zhang sv[2*k+1] = sv_wk[2*j + 1]; 4122bf73ac6SHong Zhang } 4132bf73ac6SHong Zhang } 4142bf73ac6SHong Zhang 4152bf73ac6SHong Zhang for (j=0; j<nta; j++) { 4162bf73ac6SHong Zhang ierr = PetscTableDestroy(&svtas[j]);CHKERRQ(ierr); 4172bf73ac6SHong Zhang ierr = PetscFree(ta2sv[j]);CHKERRQ(ierr); 4182bf73ac6SHong Zhang } 4192bf73ac6SHong Zhang ierr = PetscFree4(svtas,tdata,sv_wk,ta2sv);CHKERRQ(ierr); 4202bf73ac6SHong Zhang 4212bf73ac6SHong Zhang *Nsvtx = nta; 4222bf73ac6SHong Zhang *svtx = sedges; 4232bf73ac6SHong Zhang PetscFunctionReturn(0); 4242bf73ac6SHong Zhang } 4252bf73ac6SHong Zhang 426eac198afSGetnet static PetscErrorCode GetEdgelist_Coupling(DM dm,PetscInt *edges,PetscInt *nmerged_ptr,PetscInt *Nsv_ptr,SVtx **svtx_ptr) 4272bf73ac6SHong Zhang { 4282bf73ac6SHong Zhang PetscErrorCode ierr; 4292bf73ac6SHong Zhang MPI_Comm comm; 4302bf73ac6SHong Zhang PetscMPIInt size,rank,*recvcounts=NULL,*displs=NULL; 431eac198afSGetnet DM_Network *network = (DM_Network*)dm->data; 432eac198afSGetnet PetscInt i,j,ctr,np; 433eac198afSGetnet PetscInt *vidxlTog,Nsv=0,Nsubnet=network->Nsubnet; 434eac198afSGetnet PetscInt *sedgelist=network->sedgelist; 435eac198afSGetnet PetscInt net,idx,gidx,nmerged,v,*vrange,gidx_from,net_from,sv_idx; 4362bf73ac6SHong Zhang SVtxType svtype = SVNONE; 4372bf73ac6SHong Zhang SVtx *svtx=NULL; 4382bf73ac6SHong Zhang 4392bf73ac6SHong Zhang PetscFunctionBegin; 4402bf73ac6SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 44155b25c41SPierre Jolivet ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 44255b25c41SPierre Jolivet ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 4432bf73ac6SHong Zhang 4442bf73ac6SHong Zhang /* (1) Create svtx[] from sedgelist */ 4452bf73ac6SHong Zhang /* -------------------------------- */ 4462bf73ac6SHong Zhang /* Nsv: global number of SVtx; svtx: shared vertices, see SVtx in dmnetworkimpl.h */ 4472bf73ac6SHong Zhang ierr = SVtxCreate(dm,network->Nsvtx,sedgelist,&Nsv,&svtx);CHKERRQ(ierr); 4482bf73ac6SHong Zhang 4492bf73ac6SHong Zhang /* (2) Setup svtx; Shared vto vertices are merged to their vfrom vertex with same global vetex index (gidx) */ 4502bf73ac6SHong Zhang /* -------------------------------------------------------------------------------------------------------- */ 4512bf73ac6SHong Zhang /* (2.1) compute vrage[rank]: global index of 1st local vertex in proc[rank] */ 4522bf73ac6SHong Zhang ierr = PetscMalloc3(size+1,&vrange,size,&displs,size,&recvcounts);CHKERRQ(ierr); 4532bf73ac6SHong Zhang for (i=0; i<size; i++) {displs[i] = i; recvcounts[i] = 1;} 4542bf73ac6SHong Zhang 4552bf73ac6SHong Zhang vrange[0] = 0; 4562bf73ac6SHong Zhang ierr = MPI_Allgatherv(&network->nVertices,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm);CHKERRMPI(ierr); 4572bf73ac6SHong Zhang for (i=2; i<size+1; i++) { 4582bf73ac6SHong Zhang vrange[i] += vrange[i-1]; 4592bf73ac6SHong Zhang } 4602bf73ac6SHong Zhang 4612bf73ac6SHong Zhang /* (2.2) Create vidxlTog: maps UN-MERGED local vertex index i to global index gidx (plex, excluding ghost vertices) */ 4622bf73ac6SHong Zhang ierr = PetscMalloc1(network->nVertices,&vidxlTog);CHKERRQ(ierr); 4632bf73ac6SHong Zhang i = 0; gidx = 0; 4642bf73ac6SHong Zhang nmerged = 0; /* local num of merged vertices */ 4652bf73ac6SHong Zhang network->nsvtx = 0; 4662bf73ac6SHong Zhang for (net=0; net<Nsubnet; net++) { 4672bf73ac6SHong Zhang for (idx=0; idx<network->subnet[net].Nvtx; idx++) { 4682bf73ac6SHong Zhang gidx_from = gidx; 4692bf73ac6SHong Zhang sv_idx = -1; 4702bf73ac6SHong Zhang 4712bf73ac6SHong Zhang ierr = SVtxSetUp(Nsv,svtx,net,idx,&gidx_from,&svtype,&sv_idx);CHKERRQ(ierr); 4722bf73ac6SHong Zhang if (svtype == SVTO) { 4732bf73ac6SHong Zhang if (network->subnet[net].nvtx) {/* this proc owns sv_to */ 4742bf73ac6SHong Zhang net_from = svtx[sv_idx].sv[0]; /* subnet num of its shared vertex */ 4752bf73ac6SHong Zhang if (network->subnet[net_from].nvtx == 0) { 4762bf73ac6SHong Zhang /* this proc does not own v_from, thus a new local coupling vertex */ 4772bf73ac6SHong Zhang network->nsvtx++; 4782bf73ac6SHong Zhang } 4792bf73ac6SHong Zhang vidxlTog[i++] = gidx_from; 4802bf73ac6SHong Zhang nmerged++; /* a coupling vertex -- merged */ 4812bf73ac6SHong Zhang } 4822bf73ac6SHong Zhang } else { 4832bf73ac6SHong Zhang if (svtype == SVFROM) { 4842bf73ac6SHong Zhang if (network->subnet[net].nvtx) { 4852bf73ac6SHong Zhang /* this proc owns this v_from, a new local coupling vertex */ 4862bf73ac6SHong Zhang network->nsvtx++; 4872bf73ac6SHong Zhang } 4882bf73ac6SHong Zhang } 4892bf73ac6SHong Zhang if (network->subnet[net].nvtx) vidxlTog[i++] = gidx; 4902bf73ac6SHong Zhang gidx++; 4912bf73ac6SHong Zhang } 4922bf73ac6SHong Zhang } 4932bf73ac6SHong Zhang } 4942bf73ac6SHong Zhang #if defined(PETSC_USE_DEBUG) 4959ace16cdSJacob Faibussowitsch PetscAssertFalse(i != network->nVertices,PETSC_COMM_SELF,PETSC_ERR_ARG_NULL,"%D != %D nVertices",i,network->nVertices); 4962bf73ac6SHong Zhang #endif 4972bf73ac6SHong Zhang 4982bf73ac6SHong Zhang /* (2.3) Setup svtable for querry shared vertices */ 4992bf73ac6SHong Zhang for (v=0; v<Nsv; v++) { 5002bf73ac6SHong Zhang gidx = svtx[v].gidx; 5012bf73ac6SHong Zhang ierr = PetscTableAdd(network->svtable,gidx+1,v+1,INSERT_VALUES);CHKERRQ(ierr); 5022bf73ac6SHong Zhang } 5032bf73ac6SHong Zhang 5042bf73ac6SHong Zhang /* (2.4) Shared vertices in the subnetworks are merged, update global NVertices: np = sum(local nmerged) */ 5052bf73ac6SHong Zhang ierr = MPI_Allreduce(&nmerged,&np,1,MPIU_INT,MPI_SUM,comm);CHKERRMPI(ierr); 5062bf73ac6SHong Zhang network->NVertices -= np; 5072bf73ac6SHong Zhang 5082bf73ac6SHong Zhang ctr = 0; 5092bf73ac6SHong Zhang for (net=0; net<Nsubnet; net++) { 5102bf73ac6SHong Zhang for (j = 0; j < network->subnet[net].nedge; j++) { 5112bf73ac6SHong Zhang /* vfrom: */ 5122bf73ac6SHong Zhang i = network->subnet[net].edgelist[2*j] + (network->subnet[net].vStart - vrange[rank]); 5132bf73ac6SHong Zhang edges[2*ctr] = vidxlTog[i]; 5142bf73ac6SHong Zhang 5152bf73ac6SHong Zhang /* vto */ 5162bf73ac6SHong Zhang i = network->subnet[net].edgelist[2*j+1] + (network->subnet[net].vStart - vrange[rank]); 5172bf73ac6SHong Zhang edges[2*ctr+1] = vidxlTog[i]; 5182bf73ac6SHong Zhang ctr++; 5192bf73ac6SHong Zhang } 5202bf73ac6SHong Zhang } 5212bf73ac6SHong Zhang ierr = PetscFree3(vrange,displs,recvcounts);CHKERRQ(ierr); 5222bf73ac6SHong Zhang ierr = PetscFree(vidxlTog);CHKERRQ(ierr); 5232bf73ac6SHong Zhang 524eac198afSGetnet *nmerged_ptr = nmerged; 525eac198afSGetnet *Nsv_ptr = Nsv; 526eac198afSGetnet *svtx_ptr = svtx; 527eac198afSGetnet ierr = PetscFree(sedgelist);CHKERRQ(ierr); /* created in DMNetworkAddSharedVertices() */ 5285f2c45f1SShri Abhyankar PetscFunctionReturn(0); 5295f2c45f1SShri Abhyankar } 5305f2c45f1SShri Abhyankar 5315f2c45f1SShri Abhyankar /*@ 5325f2c45f1SShri Abhyankar DMNetworkLayoutSetUp - Sets up the bare layout (graph) for the network 5335f2c45f1SShri Abhyankar 534eac198afSGetnet Not Collective 5355f2c45f1SShri Abhyankar 5367a7aea1fSJed Brown Input Parameters: 537f11a936eSBarry Smith . dm - the dmnetwork object 5385f2c45f1SShri Abhyankar 5395f2c45f1SShri Abhyankar Notes: 5405f2c45f1SShri Abhyankar This routine should be called after the network sizes and edgelists have been provided. It creates 5415f2c45f1SShri Abhyankar the bare layout of the network and sets up the network to begin insertion of components. 5425f2c45f1SShri Abhyankar 5435f2c45f1SShri Abhyankar All the components should be registered before calling this routine. 5445f2c45f1SShri Abhyankar 54597bb938eSShri Abhyankar Level: beginner 5465f2c45f1SShri Abhyankar 5472bf73ac6SHong Zhang .seealso: DMNetworkSetNumSubNetworks(), DMNetworkAddSubnetwork() 5485f2c45f1SShri Abhyankar @*/ 5495f2c45f1SShri Abhyankar PetscErrorCode DMNetworkLayoutSetUp(DM dm) 5505f2c45f1SShri Abhyankar { 5515f2c45f1SShri Abhyankar PetscErrorCode ierr; 5525f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 553eac198afSGetnet PetscInt i,j,ctr,Nsubnet=network->Nsubnet,*eowners,np,*edges,*subnetvtx,*subnetedge,e,v,vfrom,vto,net; 554caf410d2SHong Zhang const PetscInt *cone; 555caf410d2SHong Zhang MPI_Comm comm; 556caf410d2SHong Zhang PetscMPIInt size,rank; 557eac198afSGetnet PetscSection sectiong; 558eac198afSGetnet PetscInt nmerged=0,Nsv=0; 559eac198afSGetnet SVtx *svtx=NULL; 5606500d4abSHong Zhang 5616500d4abSHong Zhang PetscFunctionBegin; 5629ace16cdSJacob Faibussowitsch PetscAssertFalse(network->nsubnet != network->Nsubnet,PETSC_COMM_SELF,PETSC_ERR_ARG_INCOMP,"Must call DMNetworkAddSubnetwork() %D times",network->Nsubnet); 5632bf73ac6SHong Zhang 564eac198afSGetnet /* This implementation requires user input each subnet by a single processor when Nsubnet>1, thus subnet[net].nvtx=subnet[net].Nvtx when net>0 */ 565eac198afSGetnet for (net=1; net<Nsubnet; net++) { 5669ace16cdSJacob Faibussowitsch PetscAssertFalse(network->subnet[net].nvtx && network->subnet[net].nvtx != network->subnet[net].Nvtx,PETSC_COMM_SELF,PETSC_ERR_SUP,"subnetwork %D local num of vertices %D != %D global num",net,network->subnet[net].nvtx,network->subnet[net].Nvtx); 567eac198afSGetnet } 568eac198afSGetnet 5692bf73ac6SHong Zhang /* Create svtable for querry shared vertices */ 5702bf73ac6SHong Zhang ierr = PetscTableCreate(network->Nsvtx,network->NVertices+1,&network->svtable);CHKERRQ(ierr); 5712bf73ac6SHong Zhang 572caf410d2SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 573ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 574ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 5756500d4abSHong Zhang 576f11a936eSBarry Smith /* Create LOCAL edgelist in global vertex ordering for the network by concatenating local input edgelists of the subnetworks */ 5778e71b177SVaclav Hapla ierr = PetscCalloc1(2*network->nEdges,&edges);CHKERRQ(ierr); 578eac198afSGetnet 579eac198afSGetnet if (network->Nsvtx) { /* subnetworks are coupled via shared vertices */ 580eac198afSGetnet ierr = GetEdgelist_Coupling(dm,edges,&nmerged,&Nsv,&svtx);CHKERRMPI(ierr); 581eac198afSGetnet } else { /* subnetworks are not coupled */ 582caf410d2SHong Zhang ctr = 0; 5832bf73ac6SHong Zhang for (i=0; i < Nsubnet; i++) { 5846500d4abSHong Zhang for (j = 0; j < network->subnet[i].nedge; j++) { 585ba38b151SHong Zhang edges[2*ctr] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j]; 586ba38b151SHong Zhang edges[2*ctr+1] = network->subnet[i].vStart + network->subnet[i].edgelist[2*j+1]; 5876500d4abSHong Zhang ctr++; 5886500d4abSHong Zhang } 5896500d4abSHong Zhang } 590eac198afSGetnet } 591eac198afSGetnet network->svtx = svtx; 592eac198afSGetnet network->Nsvtx = Nsv; 5937765340cSHong Zhang 5942bf73ac6SHong Zhang /* Create network->plex; One dimensional network, numCorners=2 */ 5958e71b177SVaclav Hapla ierr = DMCreate(comm,&network->plex);CHKERRQ(ierr); 5968e71b177SVaclav Hapla ierr = DMSetType(network->plex,DMPLEX);CHKERRQ(ierr); 5972bf73ac6SHong Zhang ierr = DMSetDimension(network->plex,1);CHKERRQ(ierr); 598eac198afSGetnet 599caf410d2SHong Zhang if (size == 1) { 600eac198afSGetnet ierr = DMPlexBuildFromCellList(network->plex,network->nEdges,network->nVertices-nmerged,2,edges);CHKERRQ(ierr); 601caf410d2SHong Zhang } else { 602eac198afSGetnet ierr = DMPlexBuildFromCellListParallel(network->plex,network->nEdges,network->nVertices-nmerged,PETSC_DECIDE,2,edges,NULL, NULL);CHKERRQ(ierr); 603acdb140fSBarry Smith } 6046500d4abSHong Zhang 6056500d4abSHong Zhang ierr = DMPlexGetChart(network->plex,&network->pStart,&network->pEnd);CHKERRQ(ierr); 6066500d4abSHong Zhang ierr = DMPlexGetHeightStratum(network->plex,0,&network->eStart,&network->eEnd);CHKERRQ(ierr); 6076500d4abSHong Zhang ierr = DMPlexGetHeightStratum(network->plex,1,&network->vStart,&network->vEnd);CHKERRQ(ierr); 6086500d4abSHong Zhang 609caf410d2SHong Zhang ierr = PetscSectionCreate(comm,&network->DataSection);CHKERRQ(ierr); 610caf410d2SHong Zhang ierr = PetscSectionCreate(comm,&network->DofSection);CHKERRQ(ierr); 6116500d4abSHong Zhang ierr = PetscSectionSetChart(network->DataSection,network->pStart,network->pEnd);CHKERRQ(ierr); 6126500d4abSHong Zhang ierr = PetscSectionSetChart(network->DofSection,network->pStart,network->pEnd);CHKERRQ(ierr); 6136500d4abSHong Zhang 614caf410d2SHong Zhang np = network->pEnd - network->pStart; 615caf410d2SHong Zhang ierr = PetscCalloc2(np,&network->header,np,&network->cvalue);CHKERRQ(ierr); 61654dfd506SHong Zhang for (i=0; i < np; i++) { 61754dfd506SHong Zhang network->header[i].maxcomps = 1; 61854dfd506SHong Zhang ierr = SetUpNetworkHeaderComponentValue(dm,&network->header[i],&network->cvalue[i]);CHKERRQ(ierr); 61954dfd506SHong Zhang } 620caf410d2SHong Zhang 621f11a936eSBarry Smith /* Create edge and vertex arrays for the subnetworks 622f11a936eSBarry Smith This implementation assumes that DMNetwork reads 623f11a936eSBarry Smith (1) a single subnetwork in parallel; or 624f11a936eSBarry Smith (2) n subnetworks using n processors, one subnetwork/processor. 625f11a936eSBarry Smith */ 626eac198afSGetnet ierr = PetscCalloc2(network->nEdges,&subnetedge,network->nVertices+network->nsvtx,&subnetvtx);CHKERRQ(ierr); /* Maps local edge/vertex to local subnetwork's edge/vertex */ 627f11a936eSBarry Smith network->subnetedge = subnetedge; 628f11a936eSBarry Smith network->subnetvtx = subnetvtx; 6292bf73ac6SHong Zhang for (j=0; j < network->Nsubnet; j++) { 630f11a936eSBarry Smith network->subnet[j].edges = subnetedge; 631f11a936eSBarry Smith subnetedge += network->subnet[j].nedge; 632f11a936eSBarry Smith 633f11a936eSBarry Smith network->subnet[j].vertices = subnetvtx; 634f11a936eSBarry Smith subnetvtx += network->subnet[j].nvtx; 6356500d4abSHong Zhang } 636eac198afSGetnet network->svertices = subnetvtx; 637caf410d2SHong Zhang 638caf410d2SHong Zhang /* Get edge ownership */ 6392bf73ac6SHong Zhang ierr = PetscMalloc1(size+1,&eowners);CHKERRQ(ierr); 640caf410d2SHong Zhang np = network->eEnd - network->eStart; 641ffc4695bSBarry Smith ierr = MPI_Allgather(&np,1,MPIU_INT,eowners+1,1,MPIU_INT,comm);CHKERRMPI(ierr); 642caf410d2SHong Zhang eowners[0] = 0; 643caf410d2SHong Zhang for (i=2; i<=size; i++) eowners[i] += eowners[i-1]; 644caf410d2SHong Zhang 645eac198afSGetnet /* Setup local edge and vertex arrays for subnetworks */ 6462bf73ac6SHong Zhang e = 0; 6472bf73ac6SHong Zhang for (i=0; i < Nsubnet; i++) { 648f11a936eSBarry Smith ctr = 0; 6492bf73ac6SHong Zhang for (j = 0; j < network->subnet[i].nedge; j++) { 6502bf73ac6SHong Zhang /* edge e */ 6512bf73ac6SHong Zhang network->header[e].index = e + eowners[rank]; /* Global edge index */ 6522bf73ac6SHong Zhang network->header[e].subnetid = i; 6532bf73ac6SHong Zhang network->subnet[i].edges[j] = e; 654caf410d2SHong Zhang 6552bf73ac6SHong Zhang network->header[e].ndata = 0; 6562bf73ac6SHong Zhang network->header[e].offset[0] = 0; 6572bf73ac6SHong Zhang network->header[e].offsetvarrel[0] = 0; 65854dfd506SHong Zhang ierr = PetscSectionAddDof(network->DataSection,e,network->header[e].hsize);CHKERRQ(ierr); 6592bf73ac6SHong Zhang 6602bf73ac6SHong Zhang /* connected vertices */ 6612bf73ac6SHong Zhang ierr = DMPlexGetCone(network->plex,e,&cone);CHKERRQ(ierr); 6622bf73ac6SHong Zhang 6632bf73ac6SHong Zhang /* vertex cone[0] */ 664f11a936eSBarry Smith v = cone[0]; 665f11a936eSBarry Smith network->header[v].index = edges[2*e]; /* Global vertex index */ 666f11a936eSBarry Smith network->header[v].subnetid = i; /* Subnetwork id */ 667f11a936eSBarry Smith if (Nsubnet == 1) { 668f11a936eSBarry Smith network->subnet[i].vertices[v - network->vStart] = v; /* user's subnet[].idx = petsc's v */ 669f11a936eSBarry Smith } else { 670f11a936eSBarry Smith vfrom = network->subnet[i].edgelist[2*ctr]; /* =subnet[i].idx, Global index! */ 671f11a936eSBarry Smith network->subnet[i].vertices[vfrom] = v; /* user's subnet[].dix = petsc's v */ 672f11a936eSBarry Smith } 6732bf73ac6SHong Zhang 6742bf73ac6SHong Zhang /* vertex cone[1] */ 675f11a936eSBarry Smith v = cone[1]; 676f11a936eSBarry Smith network->header[v].index = edges[2*e+1]; /* Global vertex index */ 677f11a936eSBarry Smith network->header[v].subnetid = i; /* Subnetwork id */ 678f11a936eSBarry Smith if (Nsubnet == 1) { 679f11a936eSBarry Smith network->subnet[i].vertices[v - network->vStart] = v; /* user's subnet[].idx = petsc's v */ 680f11a936eSBarry Smith } else { 681f11a936eSBarry Smith vto = network->subnet[i].edgelist[2*ctr+1]; /* =subnet[i].idx, Global index! */ 682f11a936eSBarry Smith network->subnet[i].vertices[vto] = v; /* user's subnet[].dix = petsc's v */ 683f11a936eSBarry Smith } 6842bf73ac6SHong Zhang 685f11a936eSBarry Smith e++; ctr++; 6866500d4abSHong Zhang } 6876500d4abSHong Zhang } 6882bf73ac6SHong Zhang ierr = PetscFree(eowners);CHKERRQ(ierr); 689f11a936eSBarry Smith ierr = PetscFree(edges);CHKERRQ(ierr); /* local edge list with global idx used by DMPlexBuildFromCellList() */ 690caf410d2SHong Zhang 691eac198afSGetnet /* Set local vertex array for the subnetworks */ 692eac198afSGetnet j = 0; 6932bf73ac6SHong Zhang for (v = network->vStart; v < network->vEnd; v++) { 6942bf73ac6SHong Zhang network->header[v].ndata = 0; 6952bf73ac6SHong Zhang network->header[v].offset[0] = 0; 6962bf73ac6SHong Zhang network->header[v].offsetvarrel[0] = 0; 69754dfd506SHong Zhang ierr = PetscSectionAddDof(network->DataSection,v,network->header[e].hsize);CHKERRQ(ierr); 698eac198afSGetnet 699eac198afSGetnet /* local shared vertex */ 700eac198afSGetnet ierr = PetscTableFind(network->svtable,network->header[v].index+1,&i);CHKERRQ(ierr); 701eac198afSGetnet if (i) network->svertices[j++] = v; 7026500d4abSHong Zhang } 703eac198afSGetnet 704eac198afSGetnet /* Create a global section to be used by DMNetworkIsGhostVertex() which is a non-collective routine */ 705eac198afSGetnet /* see snes_tutorials_network-ex1_4 */ 706eac198afSGetnet ierr = DMGetGlobalSection(network->plex,§iong);CHKERRQ(ierr); 7075f2c45f1SShri Abhyankar PetscFunctionReturn(0); 7085f2c45f1SShri Abhyankar } 7095f2c45f1SShri Abhyankar 71094ef8ddeSSatish Balay /*@C 7112bf73ac6SHong Zhang DMNetworkGetSubnetwork - Returns the information about a requested subnetwork 7122bf73ac6SHong Zhang 7132bf73ac6SHong Zhang Not collective 7142727e31bSShri Abhyankar 7157a7aea1fSJed Brown Input Parameters: 716caf410d2SHong Zhang + dm - the DM object 7172bf73ac6SHong Zhang - netnum - the global index of the subnetwork 7182727e31bSShri Abhyankar 7197a7aea1fSJed Brown Output Parameters: 7202727e31bSShri Abhyankar + nv - number of vertices (local) 7212727e31bSShri Abhyankar . ne - number of edges (local) 7222bf73ac6SHong Zhang . vtx - local vertices of the subnetwork 7232bf73ac6SHong Zhang - edge - local edges of the subnetwork 7242727e31bSShri Abhyankar 7252727e31bSShri Abhyankar Notes: 7262727e31bSShri Abhyankar Cannot call this routine before DMNetworkLayoutSetup() 7272727e31bSShri Abhyankar 728eac198afSGetnet The local vertices returned on each rank are determined by DMNetwork. The user does not have any control over what vertices are local. 729eac198afSGetnet 73006dd6b0eSSatish Balay Level: intermediate 73106dd6b0eSSatish Balay 7322bf73ac6SHong Zhang .seealso: DMNetworkCreate(), DMNetworkAddSubnetwork(), DMNetworkLayoutSetUp() 7332727e31bSShri Abhyankar @*/ 7342bf73ac6SHong Zhang PetscErrorCode DMNetworkGetSubnetwork(DM dm,PetscInt netnum,PetscInt *nv,PetscInt *ne,const PetscInt **vtx,const PetscInt **edge) 7352727e31bSShri Abhyankar { 736caf410d2SHong Zhang DM_Network *network = (DM_Network*)dm->data; 7372727e31bSShri Abhyankar 7382727e31bSShri Abhyankar PetscFunctionBegin; 7399ace16cdSJacob Faibussowitsch PetscAssertFalse(netnum >= network->Nsubnet,PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"Subnet index %D exceeds the num of subnets %D",netnum,network->Nsubnet); 7402bf73ac6SHong Zhang if (nv) *nv = network->subnet[netnum].nvtx; 7412bf73ac6SHong Zhang if (ne) *ne = network->subnet[netnum].nedge; 7422bf73ac6SHong Zhang if (vtx) *vtx = network->subnet[netnum].vertices; 7432bf73ac6SHong Zhang if (edge) *edge = network->subnet[netnum].edges; 7442bf73ac6SHong Zhang PetscFunctionReturn(0); 7452bf73ac6SHong Zhang } 7462bf73ac6SHong Zhang 7472bf73ac6SHong Zhang /*@ 7482bf73ac6SHong Zhang DMNetworkAddSharedVertices - Add shared vertices that connect two given subnetworks 7492bf73ac6SHong Zhang 7502bf73ac6SHong Zhang Collective on dm 7512bf73ac6SHong Zhang 7522bf73ac6SHong Zhang Input Parameters: 7532bf73ac6SHong Zhang + dm - the dm object 7542bf73ac6SHong Zhang . anetnum - first subnetwork global numbering returned by DMNetworkAddSubnetwork() 7552bf73ac6SHong Zhang . bnetnum - second subnetwork global numbering returned by DMNetworkAddSubnetwork() 7562bf73ac6SHong Zhang . nsvtx - number of vertices that are shared by the two subnetworks 7572bf73ac6SHong Zhang . asvtx - vertex index in the first subnetwork 7582bf73ac6SHong Zhang - bsvtx - vertex index in the second subnetwork 7592bf73ac6SHong Zhang 7602bf73ac6SHong Zhang Level: beginner 7612bf73ac6SHong Zhang 7622bf73ac6SHong Zhang .seealso: DMNetworkCreate(), DMNetworkAddSubnetwork(), DMNetworkGetSharedVertices() 7632bf73ac6SHong Zhang @*/ 7642bf73ac6SHong Zhang PetscErrorCode DMNetworkAddSharedVertices(DM dm,PetscInt anetnum,PetscInt bnetnum,PetscInt nsvtx,PetscInt asvtx[],PetscInt bsvtx[]) 7652bf73ac6SHong Zhang { 7662bf73ac6SHong Zhang PetscErrorCode ierr; 7672bf73ac6SHong Zhang DM_Network *network = (DM_Network*)dm->data; 7682bf73ac6SHong Zhang PetscInt i,nsubnet = network->Nsubnet,*sedgelist,Nsvtx=network->Nsvtx; 7692bf73ac6SHong Zhang 7702bf73ac6SHong Zhang PetscFunctionBegin; 7719ace16cdSJacob Faibussowitsch PetscAssertFalse(anetnum == bnetnum,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"Subnetworks must have different netnum"); 7729ace16cdSJacob Faibussowitsch PetscAssertFalse(anetnum < 0 || bnetnum < 0,PetscObjectComm((PetscObject)dm),PETSC_ERR_USER,"netnum cannot be negative"); 7732bf73ac6SHong Zhang if (!Nsvtx) { 7742bf73ac6SHong Zhang /* allocate network->sedgelist to hold at most 2*nsubnet pairs of shared vertices */ 7752bf73ac6SHong Zhang ierr = PetscMalloc1(2*4*nsubnet,&network->sedgelist);CHKERRQ(ierr); 7762bf73ac6SHong Zhang } 7772bf73ac6SHong Zhang 7782bf73ac6SHong Zhang sedgelist = network->sedgelist; 7792bf73ac6SHong Zhang for (i=0; i<nsvtx; i++) { 7802bf73ac6SHong Zhang sedgelist[4*Nsvtx] = anetnum; sedgelist[4*Nsvtx+1] = asvtx[i]; 7812bf73ac6SHong Zhang sedgelist[4*Nsvtx+2] = bnetnum; sedgelist[4*Nsvtx+3] = bsvtx[i]; 7822bf73ac6SHong Zhang Nsvtx++; 7832bf73ac6SHong Zhang } 7849ace16cdSJacob Faibussowitsch PetscAssertFalse(Nsvtx > 2*nsubnet,PETSC_COMM_SELF,PETSC_ERR_SUP,"allocate more space for coupling edgelist"); 7852bf73ac6SHong Zhang network->Nsvtx = Nsvtx; 7862727e31bSShri Abhyankar PetscFunctionReturn(0); 7872727e31bSShri Abhyankar } 7882727e31bSShri Abhyankar 7892727e31bSShri Abhyankar /*@C 7902bf73ac6SHong Zhang DMNetworkGetSharedVertices - Returns the info for the shared vertices 7912bf73ac6SHong Zhang 7922bf73ac6SHong Zhang Not collective 793caf410d2SHong Zhang 794f899ff85SJose E. Roman Input Parameter: 7952bf73ac6SHong Zhang . dm - the DM object 796caf410d2SHong Zhang 7977a7aea1fSJed Brown Output Parameters: 7982bf73ac6SHong Zhang + nsv - number of local shared vertices 7992bf73ac6SHong Zhang - svtx - local shared vertices 800caf410d2SHong Zhang 801caf410d2SHong Zhang Notes: 802caf410d2SHong Zhang Cannot call this routine before DMNetworkLayoutSetup() 803caf410d2SHong Zhang 804caf410d2SHong Zhang Level: intermediate 805caf410d2SHong Zhang 8062bf73ac6SHong Zhang .seealso: DMNetworkGetSubnetwork(), DMNetworkLayoutSetUp(), DMNetworkAddSharedVertices() 807caf410d2SHong Zhang @*/ 8082bf73ac6SHong Zhang PetscErrorCode DMNetworkGetSharedVertices(DM dm,PetscInt *nsv,const PetscInt **svtx) 809caf410d2SHong Zhang { 810caf410d2SHong Zhang DM_Network *net = (DM_Network*)dm->data; 811caf410d2SHong Zhang 812caf410d2SHong Zhang PetscFunctionBegin; 8132bf73ac6SHong Zhang if (net->Nsvtx) { 8142bf73ac6SHong Zhang *nsv = net->nsvtx; 8152bf73ac6SHong Zhang *svtx = net->svertices; 81672c3e9f7SHong Zhang } else { 8172bf73ac6SHong Zhang *nsv = 0; 8182bf73ac6SHong Zhang *svtx = NULL; 81972c3e9f7SHong Zhang } 820caf410d2SHong Zhang PetscFunctionReturn(0); 821caf410d2SHong Zhang } 822caf410d2SHong Zhang 823caf410d2SHong Zhang /*@C 8245f2c45f1SShri Abhyankar DMNetworkRegisterComponent - Registers the network component 8255f2c45f1SShri Abhyankar 826d083f849SBarry Smith Logically collective on dm 8275f2c45f1SShri Abhyankar 8287a7aea1fSJed Brown Input Parameters: 8295f2c45f1SShri Abhyankar + dm - the network object 8305f2c45f1SShri Abhyankar . name - the component name 8315f2c45f1SShri Abhyankar - size - the storage size in bytes for this component data 8325f2c45f1SShri Abhyankar 8337a7aea1fSJed Brown Output Parameters: 8345f2c45f1SShri Abhyankar . key - an integer key that defines the component 8355f2c45f1SShri Abhyankar 8365f2c45f1SShri Abhyankar Notes 8375f2c45f1SShri Abhyankar This routine should be called by all processors before calling DMNetworkLayoutSetup(). 8385f2c45f1SShri Abhyankar 83997bb938eSShri Abhyankar Level: beginner 8405f2c45f1SShri Abhyankar 8412bf73ac6SHong Zhang .seealso: DMNetworkCreate(), DMNetworkLayoutSetUp() 8425f2c45f1SShri Abhyankar @*/ 843caf410d2SHong Zhang PetscErrorCode DMNetworkRegisterComponent(DM dm,const char *name,size_t size,PetscInt *key) 8445f2c45f1SShri Abhyankar { 8455f2c45f1SShri Abhyankar PetscErrorCode ierr; 8465f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*) dm->data; 84754dfd506SHong Zhang DMNetworkComponent *component=NULL,*newcomponent=NULL; 8485f2c45f1SShri Abhyankar PetscBool flg=PETSC_FALSE; 8495f2c45f1SShri Abhyankar PetscInt i; 8505f2c45f1SShri Abhyankar 8515f2c45f1SShri Abhyankar PetscFunctionBegin; 85254dfd506SHong Zhang if (!network->component) { 85354dfd506SHong Zhang ierr = PetscCalloc1(network->max_comps_registered,&network->component);CHKERRQ(ierr); 85454dfd506SHong Zhang } 85554dfd506SHong Zhang 8565f2c45f1SShri Abhyankar for (i=0; i < network->ncomponent; i++) { 85754dfd506SHong Zhang ierr = PetscStrcmp(network->component[i].name,name,&flg);CHKERRQ(ierr); 8585f2c45f1SShri Abhyankar if (flg) { 8595f2c45f1SShri Abhyankar *key = i; 8605f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8615f2c45f1SShri Abhyankar } 8626d64e262SShri Abhyankar } 86354dfd506SHong Zhang 86454dfd506SHong Zhang if (network->ncomponent == network->max_comps_registered) { 86554dfd506SHong Zhang /* Reached max allowed so resize component */ 86654dfd506SHong Zhang network->max_comps_registered += 2; 86754dfd506SHong Zhang ierr = PetscCalloc1(network->max_comps_registered,&newcomponent);CHKERRQ(ierr); 86854dfd506SHong Zhang /* Copy over the previous component info */ 86954dfd506SHong Zhang for (i=0; i < network->ncomponent; i++) { 87054dfd506SHong Zhang ierr = PetscStrcpy(newcomponent[i].name,network->component[i].name);CHKERRQ(ierr); 87154dfd506SHong Zhang newcomponent[i].size = network->component[i].size; 8725f2c45f1SShri Abhyankar } 87354dfd506SHong Zhang /* Free old one */ 87454dfd506SHong Zhang ierr = PetscFree(network->component);CHKERRQ(ierr); 87554dfd506SHong Zhang /* Update pointer */ 87654dfd506SHong Zhang network->component = newcomponent; 87754dfd506SHong Zhang } 87854dfd506SHong Zhang 87954dfd506SHong Zhang component = &network->component[network->ncomponent]; 8805f2c45f1SShri Abhyankar 8815f2c45f1SShri Abhyankar ierr = PetscStrcpy(component->name,name);CHKERRQ(ierr); 8825f2c45f1SShri Abhyankar component->size = size/sizeof(DMNetworkComponentGenericDataType); 8835f2c45f1SShri Abhyankar *key = network->ncomponent; 8845f2c45f1SShri Abhyankar network->ncomponent++; 8855f2c45f1SShri Abhyankar PetscFunctionReturn(0); 8865f2c45f1SShri Abhyankar } 8875f2c45f1SShri Abhyankar 8885f2c45f1SShri Abhyankar /*@ 8892bf73ac6SHong Zhang DMNetworkGetVertexRange - Get the bounds [start, end) for the local vertices 8905f2c45f1SShri Abhyankar 8915f2c45f1SShri Abhyankar Not Collective 8925f2c45f1SShri Abhyankar 893f899ff85SJose E. Roman Input Parameter: 8942bf73ac6SHong Zhang . dm - the DMNetwork object 8955f2c45f1SShri Abhyankar 896fd292e60Sprj- Output Parameters: 8972bf73ac6SHong Zhang + vStart - the first vertex point 8982bf73ac6SHong Zhang - vEnd - one beyond the last vertex point 8995f2c45f1SShri Abhyankar 90097bb938eSShri Abhyankar Level: beginner 9015f2c45f1SShri Abhyankar 9022bf73ac6SHong Zhang .seealso: DMNetworkGetEdgeRange() 9035f2c45f1SShri Abhyankar @*/ 9045f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetVertexRange(DM dm,PetscInt *vStart,PetscInt *vEnd) 9055f2c45f1SShri Abhyankar { 9065f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9075f2c45f1SShri Abhyankar 9085f2c45f1SShri Abhyankar PetscFunctionBegin; 9095f2c45f1SShri Abhyankar if (vStart) *vStart = network->vStart; 9105f2c45f1SShri Abhyankar if (vEnd) *vEnd = network->vEnd; 9115f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9125f2c45f1SShri Abhyankar } 9135f2c45f1SShri Abhyankar 9145f2c45f1SShri Abhyankar /*@ 9152bf73ac6SHong Zhang DMNetworkGetEdgeRange - Get the bounds [start, end) for the local edges 9165f2c45f1SShri Abhyankar 9175f2c45f1SShri Abhyankar Not Collective 9185f2c45f1SShri Abhyankar 919f899ff85SJose E. Roman Input Parameter: 9202bf73ac6SHong Zhang . dm - the DMNetwork object 9215f2c45f1SShri Abhyankar 922fd292e60Sprj- Output Parameters: 9235f2c45f1SShri Abhyankar + eStart - The first edge point 9245f2c45f1SShri Abhyankar - eEnd - One beyond the last edge point 9255f2c45f1SShri Abhyankar 92697bb938eSShri Abhyankar Level: beginner 9275f2c45f1SShri Abhyankar 9282bf73ac6SHong Zhang .seealso: DMNetworkGetVertexRange() 9295f2c45f1SShri Abhyankar @*/ 9305f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetEdgeRange(DM dm,PetscInt *eStart,PetscInt *eEnd) 9315f2c45f1SShri Abhyankar { 9325f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 9335f2c45f1SShri Abhyankar 9345f2c45f1SShri Abhyankar PetscFunctionBegin; 9355f2c45f1SShri Abhyankar if (eStart) *eStart = network->eStart; 9365f2c45f1SShri Abhyankar if (eEnd) *eEnd = network->eEnd; 9375f2c45f1SShri Abhyankar PetscFunctionReturn(0); 9385f2c45f1SShri Abhyankar } 9395f2c45f1SShri Abhyankar 9402bf73ac6SHong Zhang static PetscErrorCode DMNetworkGetIndex(DM dm,PetscInt p,PetscInt *index) 9412bf73ac6SHong Zhang { 9422bf73ac6SHong Zhang PetscErrorCode ierr; 9432bf73ac6SHong Zhang DM_Network *network = (DM_Network*)dm->data; 9442bf73ac6SHong Zhang PetscInt offsetp; 9452bf73ac6SHong Zhang DMNetworkComponentHeader header; 9462bf73ac6SHong Zhang 9472bf73ac6SHong Zhang PetscFunctionBegin; 9489ace16cdSJacob Faibussowitsch PetscAssertFalse(!dm->setupcalled,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first"); 9492bf73ac6SHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 9502bf73ac6SHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offsetp); 9512bf73ac6SHong Zhang *index = header->index; 9522bf73ac6SHong Zhang PetscFunctionReturn(0); 9532bf73ac6SHong Zhang } 9542bf73ac6SHong Zhang 9557b6afd5bSHong Zhang /*@ 9562bf73ac6SHong Zhang DMNetworkGetGlobalEdgeIndex - Get the global numbering for the edge on the network 9577b6afd5bSHong Zhang 9587b6afd5bSHong Zhang Not Collective 9597b6afd5bSHong Zhang 9607b6afd5bSHong Zhang Input Parameters: 9617b6afd5bSHong Zhang + dm - DMNetwork object 962e85e6aecSHong Zhang - p - edge point 9637b6afd5bSHong Zhang 964fd292e60Sprj- Output Parameters: 9652bf73ac6SHong Zhang . index - the global numbering for the edge 9667b6afd5bSHong Zhang 9677b6afd5bSHong Zhang Level: intermediate 9687b6afd5bSHong Zhang 9692bf73ac6SHong Zhang .seealso: DMNetworkGetGlobalVertexIndex() 9707b6afd5bSHong Zhang @*/ 971e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalEdgeIndex(DM dm,PetscInt p,PetscInt *index) 9727b6afd5bSHong Zhang { 9737b6afd5bSHong Zhang PetscErrorCode ierr; 9747b6afd5bSHong Zhang 9757b6afd5bSHong Zhang PetscFunctionBegin; 9762bf73ac6SHong Zhang ierr = DMNetworkGetIndex(dm,p,index);CHKERRQ(ierr); 9777b6afd5bSHong Zhang PetscFunctionReturn(0); 9787b6afd5bSHong Zhang } 9797b6afd5bSHong Zhang 9805f2c45f1SShri Abhyankar /*@ 9812bf73ac6SHong Zhang DMNetworkGetGlobalVertexIndex - Get the global numbering for the vertex on the network 982e85e6aecSHong Zhang 983e85e6aecSHong Zhang Not Collective 984e85e6aecSHong Zhang 985e85e6aecSHong Zhang Input Parameters: 986e85e6aecSHong Zhang + dm - DMNetwork object 987e85e6aecSHong Zhang - p - vertex point 988e85e6aecSHong Zhang 989fd292e60Sprj- Output Parameters: 9902bf73ac6SHong Zhang . index - the global numbering for the vertex 991e85e6aecSHong Zhang 992e85e6aecSHong Zhang Level: intermediate 993e85e6aecSHong Zhang 9942bf73ac6SHong Zhang .seealso: DMNetworkGetGlobalEdgeIndex() 995e85e6aecSHong Zhang @*/ 996e85e6aecSHong Zhang PetscErrorCode DMNetworkGetGlobalVertexIndex(DM dm,PetscInt p,PetscInt *index) 997e85e6aecSHong Zhang { 998e85e6aecSHong Zhang PetscErrorCode ierr; 999e85e6aecSHong Zhang 1000e85e6aecSHong Zhang PetscFunctionBegin; 10012bf73ac6SHong Zhang ierr = DMNetworkGetIndex(dm,p,index);CHKERRQ(ierr); 10029988b915SShri Abhyankar PetscFunctionReturn(0); 10039988b915SShri Abhyankar } 10049988b915SShri Abhyankar 10059988b915SShri Abhyankar /*@ 10065f2c45f1SShri Abhyankar DMNetworkGetNumComponents - Get the number of components at a vertex/edge 10075f2c45f1SShri Abhyankar 10085f2c45f1SShri Abhyankar Not Collective 10095f2c45f1SShri Abhyankar 10105f2c45f1SShri Abhyankar Input Parameters: 10112bf73ac6SHong Zhang + dm - the DMNetwork object 1012a2b725a8SWilliam Gropp - p - vertex/edge point 10135f2c45f1SShri Abhyankar 10145f2c45f1SShri Abhyankar Output Parameters: 10155f2c45f1SShri Abhyankar . numcomponents - Number of components at the vertex/edge 10165f2c45f1SShri Abhyankar 101797bb938eSShri Abhyankar Level: beginner 10185f2c45f1SShri Abhyankar 10192bf73ac6SHong Zhang .seealso: DMNetworkRegisterComponent(), DMNetworkAddComponent() 10205f2c45f1SShri Abhyankar @*/ 10215f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetNumComponents(DM dm,PetscInt p,PetscInt *numcomponents) 10225f2c45f1SShri Abhyankar { 10235f2c45f1SShri Abhyankar PetscErrorCode ierr; 10245f2c45f1SShri Abhyankar PetscInt offset; 10255f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 10265f2c45f1SShri Abhyankar 10275f2c45f1SShri Abhyankar PetscFunctionBegin; 10285f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 10295f2c45f1SShri Abhyankar *numcomponents = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 10305f2c45f1SShri Abhyankar PetscFunctionReturn(0); 10315f2c45f1SShri Abhyankar } 10325f2c45f1SShri Abhyankar 10335f2c45f1SShri Abhyankar /*@ 10342bf73ac6SHong Zhang DMNetworkGetLocalVecOffset - Get the offset for accessing the variables associated with a component at the given vertex/edge from the local vector 10355f2c45f1SShri Abhyankar 10365f2c45f1SShri Abhyankar Not Collective 10375f2c45f1SShri Abhyankar 10385f2c45f1SShri Abhyankar Input Parameters: 10392bf73ac6SHong Zhang + dm - the DMNetwork object 10407d928bffSSatish Balay . p - the edge/vertex point 10412bf73ac6SHong Zhang - compnum - component number; use ALL_COMPONENTS if no specific component is requested 10429988b915SShri Abhyankar 10439988b915SShri Abhyankar Output Parameters: 10442bf73ac6SHong Zhang . offset - the local offset 10459988b915SShri Abhyankar 10469988b915SShri Abhyankar Level: intermediate 10479988b915SShri Abhyankar 10482bf73ac6SHong Zhang .seealso: DMGetLocalVector(), DMNetworkGetComponent(), DMNetworkGetGlobalVecOffset() 10499988b915SShri Abhyankar @*/ 10502bf73ac6SHong Zhang PetscErrorCode DMNetworkGetLocalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offset) 10519988b915SShri Abhyankar { 10529988b915SShri Abhyankar PetscErrorCode ierr; 10539988b915SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 10549988b915SShri Abhyankar PetscInt offsetp,offsetd; 10559988b915SShri Abhyankar DMNetworkComponentHeader header; 10569988b915SShri Abhyankar 10579988b915SShri Abhyankar PetscFunctionBegin; 10582bf73ac6SHong Zhang ierr = PetscSectionGetOffset(network->plex->localSection,p,&offsetp);CHKERRQ(ierr); 10592bf73ac6SHong Zhang if (compnum == ALL_COMPONENTS) { 10602bf73ac6SHong Zhang *offset = offsetp; 10612bf73ac6SHong Zhang PetscFunctionReturn(0); 10622bf73ac6SHong Zhang } 10632bf73ac6SHong Zhang 10649988b915SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetd);CHKERRQ(ierr); 10659988b915SShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd); 10669988b915SShri Abhyankar *offset = offsetp + header->offsetvarrel[compnum]; 10679988b915SShri Abhyankar PetscFunctionReturn(0); 10689988b915SShri Abhyankar } 10699988b915SShri Abhyankar 10709988b915SShri Abhyankar /*@ 10712bf73ac6SHong Zhang DMNetworkGetGlobalVecOffset - Get the global offset for accessing the variables associated with a component for the given vertex/edge from the global vector 10729988b915SShri Abhyankar 10739988b915SShri Abhyankar Not Collective 10749988b915SShri Abhyankar 10759988b915SShri Abhyankar Input Parameters: 10762bf73ac6SHong Zhang + dm - the DMNetwork object 10777d928bffSSatish Balay . p - the edge/vertex point 10782bf73ac6SHong Zhang - compnum - component number; use ALL_COMPONENTS if no specific component is requested 10799988b915SShri Abhyankar 10809988b915SShri Abhyankar Output Parameters: 10819988b915SShri Abhyankar . offsetg - the global offset 10829988b915SShri Abhyankar 10839988b915SShri Abhyankar Level: intermediate 10849988b915SShri Abhyankar 10852bf73ac6SHong Zhang .seealso: DMNetworkGetLocalVecOffset(), DMGetGlobalVector(), DMNetworkGetComponent() 10869988b915SShri Abhyankar @*/ 10872bf73ac6SHong Zhang PetscErrorCode DMNetworkGetGlobalVecOffset(DM dm,PetscInt p,PetscInt compnum,PetscInt *offsetg) 10889988b915SShri Abhyankar { 10899988b915SShri Abhyankar PetscErrorCode ierr; 10909988b915SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 10919988b915SShri Abhyankar PetscInt offsetp,offsetd; 10929988b915SShri Abhyankar DMNetworkComponentHeader header; 10939988b915SShri Abhyankar 10949988b915SShri Abhyankar PetscFunctionBegin; 10952bf73ac6SHong Zhang ierr = PetscSectionGetOffset(network->plex->globalSection,p,&offsetp);CHKERRQ(ierr); 10962bf73ac6SHong Zhang if (offsetp < 0) offsetp = -(offsetp + 1); /* Convert to actual global offset for ghost vertex */ 10972bf73ac6SHong Zhang 10982bf73ac6SHong Zhang if (compnum == ALL_COMPONENTS) { 10992bf73ac6SHong Zhang *offsetg = offsetp; 11002bf73ac6SHong Zhang PetscFunctionReturn(0); 11012bf73ac6SHong Zhang } 11029988b915SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetd);CHKERRQ(ierr); 11039988b915SShri Abhyankar header = (DMNetworkComponentHeader)(network->componentdataarray+offsetd); 11049988b915SShri Abhyankar *offsetg = offsetp + header->offsetvarrel[compnum]; 11059988b915SShri Abhyankar PetscFunctionReturn(0); 11069988b915SShri Abhyankar } 11079988b915SShri Abhyankar 11089988b915SShri Abhyankar /*@ 11092bf73ac6SHong Zhang DMNetworkGetEdgeOffset - Get the offset for accessing the variables associated with the given edge from the local subvector 111024121865SAdrian Maldonado 111124121865SAdrian Maldonado Not Collective 111224121865SAdrian Maldonado 111324121865SAdrian Maldonado Input Parameters: 11142bf73ac6SHong Zhang + dm - the DMNetwork object 111524121865SAdrian Maldonado - p - the edge point 111624121865SAdrian Maldonado 111724121865SAdrian Maldonado Output Parameters: 111824121865SAdrian Maldonado . offset - the offset 111924121865SAdrian Maldonado 112024121865SAdrian Maldonado Level: intermediate 112124121865SAdrian Maldonado 11222bf73ac6SHong Zhang .seealso: DMNetworkGetLocalVecOffset(), DMGetLocalVector() 112324121865SAdrian Maldonado @*/ 112424121865SAdrian Maldonado PetscErrorCode DMNetworkGetEdgeOffset(DM dm,PetscInt p,PetscInt *offset) 112524121865SAdrian Maldonado { 112624121865SAdrian Maldonado PetscErrorCode ierr; 112724121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 112824121865SAdrian Maldonado 112924121865SAdrian Maldonado PetscFunctionBegin; 113024121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->edge.DofSection,p,offset);CHKERRQ(ierr); 113124121865SAdrian Maldonado PetscFunctionReturn(0); 113224121865SAdrian Maldonado } 113324121865SAdrian Maldonado 113424121865SAdrian Maldonado /*@ 11352bf73ac6SHong Zhang DMNetworkGetVertexOffset - Get the offset for accessing the variables associated with the given vertex from the local subvector 113624121865SAdrian Maldonado 113724121865SAdrian Maldonado Not Collective 113824121865SAdrian Maldonado 113924121865SAdrian Maldonado Input Parameters: 11402bf73ac6SHong Zhang + dm - the DMNetwork object 114124121865SAdrian Maldonado - p - the vertex point 114224121865SAdrian Maldonado 114324121865SAdrian Maldonado Output Parameters: 114424121865SAdrian Maldonado . offset - the offset 114524121865SAdrian Maldonado 114624121865SAdrian Maldonado Level: intermediate 114724121865SAdrian Maldonado 11482bf73ac6SHong Zhang .seealso: DMNetworkGetEdgeOffset(), DMGetLocalVector() 114924121865SAdrian Maldonado @*/ 115024121865SAdrian Maldonado PetscErrorCode DMNetworkGetVertexOffset(DM dm,PetscInt p,PetscInt *offset) 115124121865SAdrian Maldonado { 115224121865SAdrian Maldonado PetscErrorCode ierr; 115324121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 115424121865SAdrian Maldonado 115524121865SAdrian Maldonado PetscFunctionBegin; 115624121865SAdrian Maldonado p -= network->vStart; 115724121865SAdrian Maldonado ierr = PetscSectionGetOffset(network->vertex.DofSection,p,offset);CHKERRQ(ierr); 115824121865SAdrian Maldonado PetscFunctionReturn(0); 115924121865SAdrian Maldonado } 11602bf73ac6SHong Zhang 11615f2c45f1SShri Abhyankar /*@ 11622bf73ac6SHong Zhang DMNetworkAddComponent - Adds a network component and number of variables at the given point (vertex/edge) 11635f2c45f1SShri Abhyankar 1164eac198afSGetnet Collective on dm 11655f2c45f1SShri Abhyankar 11665f2c45f1SShri Abhyankar Input Parameters: 11674dc646bcSBarry Smith + dm - the DMNetwork 1168eac198afSGetnet . p - the vertex/edge point. These points are local indices provided by DMNetworkGetSubnetwork() 1169eac198afSGetnet . componentkey - component key returned while registering the component with DMNetworkRegisterComponent() 1170eac198afSGetnet . compvalue - pointer to the data structure for the component, or NULL if the component does not require data 1171eac198afSGetnet - nvar - number of variables for the component at the vertex/edge point, zero if the component does not introduce any degrees of freedom at the point 11725f2c45f1SShri Abhyankar 1173eac198afSGetnet Notes: 1174eac198afSGetnet The owning rank and any other ranks that have this point as a ghost location must call this routine to add a component and number of variables in the same order at the given point. 1175eac198afSGetnet 1176eac198afSGetnet DMNetworkLayoutSetUp() must be called before this routine. 1177eac198afSGetnet 1178eac198afSGetnet Developer Notes: 1179eac198afSGetnet The requirement that all the ranks with access to a vertex (as owner or as ghost) add all the components comes from a limitation of the underlying implementation based on DMPLEX. 118097bb938eSShri Abhyankar Level: beginner 11815f2c45f1SShri Abhyankar 1182eac198afSGetnet .seealso: DMNetworkGetComponent(), DMNetworkGetSubnetwork(), DMNetworkIsGhostVertex(), DMNetworkLayoutSetUp() 11835f2c45f1SShri Abhyankar @*/ 11842bf73ac6SHong Zhang PetscErrorCode DMNetworkAddComponent(DM dm,PetscInt p,PetscInt componentkey,void* compvalue,PetscInt nvar) 11855f2c45f1SShri Abhyankar { 11865f2c45f1SShri Abhyankar PetscErrorCode ierr; 11875f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 11882bf73ac6SHong Zhang DMNetworkComponent *component = &network->component[componentkey]; 118954dfd506SHong Zhang DMNetworkComponentHeader header; 119054dfd506SHong Zhang DMNetworkComponentValue cvalue; 119154dfd506SHong Zhang PetscInt compnum; 119254dfd506SHong Zhang PetscInt *compsize,*compkey,*compoffset,*compnvar,*compoffsetvarrel; 119354dfd506SHong Zhang void* *compdata; 11945f2c45f1SShri Abhyankar 11955f2c45f1SShri Abhyankar PetscFunctionBegin; 11969ace16cdSJacob Faibussowitsch PetscAssertFalse(componentkey < 0,PETSC_COMM_SELF,PETSC_ERR_ARG_WRONG,"componentkey %D cannot be negative. Input a component key returned while registering the component with DMNetworkRegisterComponent()",componentkey); 1197eac198afSGetnet 1198eac198afSGetnet /* The owning rank and all ghost ranks add nvar */ 11995f2c45f1SShri Abhyankar ierr = PetscSectionAddDof(network->DofSection,p,nvar);CHKERRQ(ierr); 12002bf73ac6SHong Zhang 1201eac198afSGetnet /* The owning rank and all ghost ranks add a component, including compvalue=NULL */ 120254dfd506SHong Zhang header = &network->header[p]; 120354dfd506SHong Zhang cvalue = &network->cvalue[p]; 120454dfd506SHong Zhang if (header->ndata == header->maxcomps) { 1205f11a936eSBarry Smith PetscInt additional_size; 1206f11a936eSBarry Smith 120754dfd506SHong Zhang /* Reached limit so resize header component arrays */ 120854dfd506SHong Zhang header->maxcomps += 2; 120954dfd506SHong Zhang 121054dfd506SHong Zhang /* Allocate arrays for component information and value */ 121154dfd506SHong Zhang ierr = PetscCalloc5(header->maxcomps,&compsize,header->maxcomps,&compkey,header->maxcomps,&compoffset,header->maxcomps,&compnvar,header->maxcomps,&compoffsetvarrel);CHKERRQ(ierr); 121254dfd506SHong Zhang ierr = PetscMalloc1(header->maxcomps,&compdata);CHKERRQ(ierr); 121354dfd506SHong Zhang 121454dfd506SHong Zhang /* Recalculate header size */ 121554dfd506SHong Zhang header->hsize = sizeof(struct _p_DMNetworkComponentHeader) + 5*header->maxcomps*sizeof(PetscInt); 121654dfd506SHong Zhang 121754dfd506SHong Zhang header->hsize /= sizeof(DMNetworkComponentGenericDataType); 121854dfd506SHong Zhang 121954dfd506SHong Zhang /* Copy over component info */ 122054dfd506SHong Zhang ierr = PetscMemcpy(compsize,header->size,header->ndata*sizeof(PetscInt));CHKERRQ(ierr); 122154dfd506SHong Zhang ierr = PetscMemcpy(compkey,header->key,header->ndata*sizeof(PetscInt));CHKERRQ(ierr); 122254dfd506SHong Zhang ierr = PetscMemcpy(compoffset,header->offset,header->ndata*sizeof(PetscInt));CHKERRQ(ierr); 122354dfd506SHong Zhang ierr = PetscMemcpy(compnvar,header->nvar,header->ndata*sizeof(PetscInt));CHKERRQ(ierr); 122454dfd506SHong Zhang ierr = PetscMemcpy(compoffsetvarrel,header->offsetvarrel,header->ndata*sizeof(PetscInt));CHKERRQ(ierr); 122554dfd506SHong Zhang 122654dfd506SHong Zhang /* Copy over component data pointers */ 122754dfd506SHong Zhang ierr = PetscMemcpy(compdata,cvalue->data,header->ndata*sizeof(void*));CHKERRQ(ierr); 122854dfd506SHong Zhang 122954dfd506SHong Zhang /* Free old arrays */ 123054dfd506SHong Zhang ierr = PetscFree5(header->size,header->key,header->offset,header->nvar,header->offsetvarrel);CHKERRQ(ierr); 123154dfd506SHong Zhang ierr = PetscFree(cvalue->data);CHKERRQ(ierr); 123254dfd506SHong Zhang 123354dfd506SHong Zhang /* Update pointers */ 123454dfd506SHong Zhang header->size = compsize; 123554dfd506SHong Zhang header->key = compkey; 123654dfd506SHong Zhang header->offset = compoffset; 123754dfd506SHong Zhang header->nvar = compnvar; 123854dfd506SHong Zhang header->offsetvarrel = compoffsetvarrel; 123954dfd506SHong Zhang 124054dfd506SHong Zhang cvalue->data = compdata; 124154dfd506SHong Zhang 124254dfd506SHong Zhang /* Update DataSection Dofs */ 124354dfd506SHong Zhang /* The dofs for datasection point p equals sizeof the header (i.e. header->hsize) + sizes of the components added at point p. With the resizing of the header, we need to update the dofs for point p. Hence, we add the extra size added for the header */ 1244f11a936eSBarry Smith additional_size = (5*(header->maxcomps - header->ndata)*sizeof(PetscInt))/sizeof(DMNetworkComponentGenericDataType); 124554dfd506SHong Zhang ierr = PetscSectionAddDof(network->DataSection,p,additional_size);CHKERRQ(ierr); 124654dfd506SHong Zhang } 124754dfd506SHong Zhang header = &network->header[p]; 124854dfd506SHong Zhang cvalue = &network->cvalue[p]; 124954dfd506SHong Zhang 125054dfd506SHong Zhang compnum = header->ndata; 12512bf73ac6SHong Zhang 12522bf73ac6SHong Zhang header->size[compnum] = component->size; 12532bf73ac6SHong Zhang ierr = PetscSectionAddDof(network->DataSection,p,component->size);CHKERRQ(ierr); 12542bf73ac6SHong Zhang header->key[compnum] = componentkey; 12552bf73ac6SHong Zhang if (compnum != 0) header->offset[compnum] = header->offset[compnum-1] + header->size[compnum-1]; 12562bf73ac6SHong Zhang cvalue->data[compnum] = (void*)compvalue; 12572bf73ac6SHong Zhang 12582bf73ac6SHong Zhang /* variables */ 12592bf73ac6SHong Zhang header->nvar[compnum] += nvar; 12602bf73ac6SHong Zhang if (compnum != 0) header->offsetvarrel[compnum] = header->offsetvarrel[compnum-1] + header->nvar[compnum-1]; 12612bf73ac6SHong Zhang 12622bf73ac6SHong Zhang header->ndata++; 12635f2c45f1SShri Abhyankar PetscFunctionReturn(0); 12645f2c45f1SShri Abhyankar } 12655f2c45f1SShri Abhyankar 126627f51fceSHong Zhang /*@ 12672bf73ac6SHong Zhang DMNetworkGetComponent - Gets the component key, the component data, and the number of variables at a given network point 126827f51fceSHong Zhang 126927f51fceSHong Zhang Not Collective 127027f51fceSHong Zhang 127127f51fceSHong Zhang Input Parameters: 12722bf73ac6SHong Zhang + dm - the DMNetwork object 12732bf73ac6SHong Zhang . p - vertex/edge point 127499c90e12SSatish Balay - compnum - component number; use ALL_COMPONENTS if sum up all the components 127527f51fceSHong Zhang 127627f51fceSHong Zhang Output Parameters: 12772bf73ac6SHong Zhang + compkey - the key obtained when registering the component (use NULL if not required) 12782bf73ac6SHong Zhang . component - the component data (use NULL if not required) 12792bf73ac6SHong Zhang - nvar - number of variables (use NULL if not required) 128027f51fceSHong Zhang 128197bb938eSShri Abhyankar Level: beginner 128227f51fceSHong Zhang 12832bf73ac6SHong Zhang .seealso: DMNetworkAddComponent(), DMNetworkGetNumComponents() 128427f51fceSHong Zhang @*/ 12852bf73ac6SHong Zhang PetscErrorCode DMNetworkGetComponent(DM dm,PetscInt p,PetscInt compnum,PetscInt *compkey,void **component,PetscInt *nvar) 128627f51fceSHong Zhang { 128727f51fceSHong Zhang PetscErrorCode ierr; 128827f51fceSHong Zhang DM_Network *network = (DM_Network*)dm->data; 12892bf73ac6SHong Zhang PetscInt offset = 0; 12902bf73ac6SHong Zhang DMNetworkComponentHeader header; 129127f51fceSHong Zhang 129227f51fceSHong Zhang PetscFunctionBegin; 12932bf73ac6SHong Zhang if (compnum == ALL_COMPONENTS) { 129427f51fceSHong Zhang ierr = PetscSectionGetDof(network->DofSection,p,nvar);CHKERRQ(ierr); 129527f51fceSHong Zhang PetscFunctionReturn(0); 129627f51fceSHong Zhang } 129727f51fceSHong Zhang 12982bf73ac6SHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 129942dc13f1SHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offset); 13005f2c45f1SShri Abhyankar 13012bf73ac6SHong Zhang if (compnum >= 0) { 13022bf73ac6SHong Zhang if (compkey) *compkey = header->key[compnum]; 13032bf73ac6SHong Zhang if (component) { 130454dfd506SHong Zhang offset += header->hsize+header->offset[compnum]; 13052bf73ac6SHong Zhang *component = network->componentdataarray+offset; 13062bf73ac6SHong Zhang } 13072bf73ac6SHong Zhang } 13085f2c45f1SShri Abhyankar 13092bf73ac6SHong Zhang if (nvar) *nvar = header->nvar[compnum]; 131054dfd506SHong Zhang 13115f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13125f2c45f1SShri Abhyankar } 13135f2c45f1SShri Abhyankar 13142bf73ac6SHong Zhang /* 13152bf73ac6SHong Zhang Sets up the array that holds the data for all components and its associated section. 13162bf73ac6SHong Zhang It copies the data for all components in a contiguous array called componentdataarray. The component data is stored pointwise with an additional header (metadata) stored for each point. The header has metadata information such as number of components at each point, number of variables for each component, offsets for the components data, etc. 13172bf73ac6SHong Zhang */ 13185f2c45f1SShri Abhyankar PetscErrorCode DMNetworkComponentSetUp(DM dm) 13195f2c45f1SShri Abhyankar { 13205f2c45f1SShri Abhyankar PetscErrorCode ierr; 13215f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 1322f11a936eSBarry Smith PetscInt arr_size,p,offset,offsetp,ncomp,i,*headerarr; 13235f2c45f1SShri Abhyankar DMNetworkComponentHeader header; 13245f2c45f1SShri Abhyankar DMNetworkComponentValue cvalue; 1325f11a936eSBarry Smith DMNetworkComponentHeader headerinfo; 13265f2c45f1SShri Abhyankar DMNetworkComponentGenericDataType *componentdataarray; 13275f2c45f1SShri Abhyankar 13285f2c45f1SShri Abhyankar PetscFunctionBegin; 13295f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DataSection);CHKERRQ(ierr); 13305f2c45f1SShri Abhyankar ierr = PetscSectionGetStorageSize(network->DataSection,&arr_size);CHKERRQ(ierr); 1331f11a936eSBarry Smith /* arr_size+1 fixes pipeline test of opensolaris-misc for src/dm/tests/ex10.c -- Do not know why */ 1332f11a936eSBarry Smith ierr = PetscCalloc1(arr_size+1,&network->componentdataarray);CHKERRQ(ierr); 1333eac198afSGetnet 13345f2c45f1SShri Abhyankar componentdataarray = network->componentdataarray; 13355f2c45f1SShri Abhyankar for (p = network->pStart; p < network->pEnd; p++) { 13365f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(network->DataSection,p,&offsetp);CHKERRQ(ierr); 13375f2c45f1SShri Abhyankar /* Copy header */ 13385f2c45f1SShri Abhyankar header = &network->header[p]; 1339f11a936eSBarry Smith headerinfo = (DMNetworkComponentHeader)(componentdataarray+offsetp); 134054dfd506SHong Zhang ierr = PetscMemcpy(headerinfo,header,sizeof(struct _p_DMNetworkComponentHeader));CHKERRQ(ierr); 1341f11a936eSBarry Smith headerarr = (PetscInt*)(headerinfo+1); 134254dfd506SHong Zhang ierr = PetscMemcpy(headerarr,header->size,header->maxcomps*sizeof(PetscInt));CHKERRQ(ierr); 134354dfd506SHong Zhang headerarr += header->maxcomps; 134454dfd506SHong Zhang ierr = PetscMemcpy(headerarr,header->key,header->maxcomps*sizeof(PetscInt));CHKERRQ(ierr); 134554dfd506SHong Zhang headerarr += header->maxcomps; 134654dfd506SHong Zhang ierr = PetscMemcpy(headerarr,header->offset,header->maxcomps*sizeof(PetscInt));CHKERRQ(ierr); 134754dfd506SHong Zhang headerarr += header->maxcomps; 134854dfd506SHong Zhang ierr = PetscMemcpy(headerarr,header->nvar,header->maxcomps*sizeof(PetscInt));CHKERRQ(ierr); 134954dfd506SHong Zhang headerarr += header->maxcomps; 135054dfd506SHong Zhang ierr = PetscMemcpy(headerarr,header->offsetvarrel,header->maxcomps*sizeof(PetscInt));CHKERRQ(ierr); 135154dfd506SHong Zhang 13525f2c45f1SShri Abhyankar /* Copy data */ 13535f2c45f1SShri Abhyankar cvalue = &network->cvalue[p]; 13545f2c45f1SShri Abhyankar ncomp = header->ndata; 13552bf73ac6SHong Zhang 13565f2c45f1SShri Abhyankar for (i = 0; i < ncomp; i++) { 135754dfd506SHong Zhang offset = offsetp + header->hsize + header->offset[i]; 1358302440fdSBarry Smith ierr = PetscMemcpy(componentdataarray+offset,cvalue->data[i],header->size[i]*sizeof(DMNetworkComponentGenericDataType));CHKERRQ(ierr); 13595f2c45f1SShri Abhyankar } 13605f2c45f1SShri Abhyankar } 13615f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13625f2c45f1SShri Abhyankar } 13635f2c45f1SShri Abhyankar 13645f2c45f1SShri Abhyankar /* Sets up the section for dofs. This routine is called during DMSetUp() */ 13652bf73ac6SHong Zhang static PetscErrorCode DMNetworkVariablesSetUp(DM dm) 13665f2c45f1SShri Abhyankar { 13675f2c45f1SShri Abhyankar PetscErrorCode ierr; 13685f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 13695f2c45f1SShri Abhyankar 13705f2c45f1SShri Abhyankar PetscFunctionBegin; 13715f2c45f1SShri Abhyankar ierr = PetscSectionSetUp(network->DofSection);CHKERRQ(ierr); 13725f2c45f1SShri Abhyankar PetscFunctionReturn(0); 13735f2c45f1SShri Abhyankar } 13745f2c45f1SShri Abhyankar 137524121865SAdrian Maldonado /* Get a subsection from a range of points */ 13769dddd249SSatish Balay static PetscErrorCode DMNetworkGetSubSection_private(PetscSection main,PetscInt pstart,PetscInt pend,PetscSection *subsection) 137724121865SAdrian Maldonado { 137824121865SAdrian Maldonado PetscErrorCode ierr; 137924121865SAdrian Maldonado PetscInt i, nvar; 138024121865SAdrian Maldonado 138124121865SAdrian Maldonado PetscFunctionBegin; 13829dddd249SSatish Balay ierr = PetscSectionCreate(PetscObjectComm((PetscObject)main), subsection);CHKERRQ(ierr); 138324121865SAdrian Maldonado ierr = PetscSectionSetChart(*subsection, 0, pend - pstart);CHKERRQ(ierr); 138424121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 13859dddd249SSatish Balay ierr = PetscSectionGetDof(main,i,&nvar);CHKERRQ(ierr); 138624121865SAdrian Maldonado ierr = PetscSectionSetDof(*subsection, i - pstart, nvar);CHKERRQ(ierr); 138724121865SAdrian Maldonado } 138824121865SAdrian Maldonado 138924121865SAdrian Maldonado ierr = PetscSectionSetUp(*subsection);CHKERRQ(ierr); 139024121865SAdrian Maldonado PetscFunctionReturn(0); 139124121865SAdrian Maldonado } 139224121865SAdrian Maldonado 139324121865SAdrian Maldonado /* Create a submap of points with a GlobalToLocal structure */ 13942bf73ac6SHong Zhang static PetscErrorCode DMNetworkSetSubMap_private(PetscInt pstart, PetscInt pend, ISLocalToGlobalMapping *map) 139524121865SAdrian Maldonado { 139624121865SAdrian Maldonado PetscErrorCode ierr; 139724121865SAdrian Maldonado PetscInt i, *subpoints; 139824121865SAdrian Maldonado 139924121865SAdrian Maldonado PetscFunctionBegin; 140024121865SAdrian Maldonado /* Create index sets to map from "points" to "subpoints" */ 140124121865SAdrian Maldonado ierr = PetscMalloc1(pend - pstart, &subpoints);CHKERRQ(ierr); 140224121865SAdrian Maldonado for (i = pstart; i < pend; i++) { 140324121865SAdrian Maldonado subpoints[i - pstart] = i; 140424121865SAdrian Maldonado } 1405459726d8SSatish Balay ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,pend-pstart,subpoints,PETSC_COPY_VALUES,map);CHKERRQ(ierr); 140624121865SAdrian Maldonado ierr = PetscFree(subpoints);CHKERRQ(ierr); 140724121865SAdrian Maldonado PetscFunctionReturn(0); 140824121865SAdrian Maldonado } 140924121865SAdrian Maldonado 141024121865SAdrian Maldonado /*@ 14112bf73ac6SHong Zhang DMNetworkAssembleGraphStructures - Assembles vertex and edge data structures. Must be called after DMNetworkDistribute 141224121865SAdrian Maldonado 14132bf73ac6SHong Zhang Collective on dm 141424121865SAdrian Maldonado 141524121865SAdrian Maldonado Input Parameters: 14162bf73ac6SHong Zhang . dm - the DMNetworkObject 141724121865SAdrian Maldonado 141824121865SAdrian Maldonado Note: the routine will create alternative orderings for the vertices and edges. Assume global network points are: 141924121865SAdrian Maldonado 142024121865SAdrian Maldonado points = [0 1 2 3 4 5 6] 142124121865SAdrian Maldonado 1422caf410d2SHong 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]). 142324121865SAdrian Maldonado 142424121865SAdrian Maldonado With this new ordering a local PetscSection, global PetscSection and PetscSF will be created specific to the subset. 142524121865SAdrian Maldonado 142624121865SAdrian Maldonado Level: intermediate 142724121865SAdrian Maldonado 142824121865SAdrian Maldonado @*/ 142924121865SAdrian Maldonado PetscErrorCode DMNetworkAssembleGraphStructures(DM dm) 143024121865SAdrian Maldonado { 143124121865SAdrian Maldonado PetscErrorCode ierr; 143224121865SAdrian Maldonado MPI_Comm comm; 1433f11a936eSBarry Smith PetscMPIInt size; 143424121865SAdrian Maldonado DM_Network *network = (DM_Network*)dm->data; 143524121865SAdrian Maldonado 1436eab1376dSHong Zhang PetscFunctionBegin; 143724121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 1438ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 143924121865SAdrian Maldonado 144024121865SAdrian Maldonado /* Create maps for vertices and edges */ 144124121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 144224121865SAdrian Maldonado ierr = DMNetworkSetSubMap_private(network->eStart,network->eEnd,&network->edge.mapping);CHKERRQ(ierr); 144324121865SAdrian Maldonado 144424121865SAdrian Maldonado /* Create local sub-sections */ 144524121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->vStart,network->vEnd,&network->vertex.DofSection);CHKERRQ(ierr); 144624121865SAdrian Maldonado ierr = DMNetworkGetSubSection_private(network->DofSection,network->eStart,network->eEnd,&network->edge.DofSection);CHKERRQ(ierr); 144724121865SAdrian Maldonado 14489852e123SBarry Smith if (size > 1) { 144924121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 145022bbedd7SHong Zhang 145124121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->vertex.DofSection, network->vertex.sf, PETSC_FALSE, PETSC_FALSE, &network->vertex.GlobalDofSection);CHKERRQ(ierr); 145224121865SAdrian Maldonado ierr = PetscSFGetSubSF(network->plex->sf, network->edge.mapping, &network->edge.sf);CHKERRQ(ierr); 145324121865SAdrian Maldonado ierr = PetscSectionCreateGlobalSection(network->edge.DofSection, network->edge.sf, PETSC_FALSE, PETSC_FALSE, &network->edge.GlobalDofSection);CHKERRQ(ierr); 145424121865SAdrian Maldonado } else { 145524121865SAdrian Maldonado /* create structures for vertex */ 145624121865SAdrian Maldonado ierr = PetscSectionClone(network->vertex.DofSection,&network->vertex.GlobalDofSection);CHKERRQ(ierr); 145724121865SAdrian Maldonado /* create structures for edge */ 145824121865SAdrian Maldonado ierr = PetscSectionClone(network->edge.DofSection,&network->edge.GlobalDofSection);CHKERRQ(ierr); 145924121865SAdrian Maldonado } 146024121865SAdrian Maldonado 146124121865SAdrian Maldonado /* Add viewers */ 146224121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->edge.GlobalDofSection,"Global edge dof section");CHKERRQ(ierr); 146324121865SAdrian Maldonado ierr = PetscObjectSetName((PetscObject)network->vertex.GlobalDofSection,"Global vertex dof section");CHKERRQ(ierr); 146424121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->edge.GlobalDofSection, NULL, "-edge_global_section_view");CHKERRQ(ierr); 146524121865SAdrian Maldonado ierr = PetscSectionViewFromOptions(network->vertex.GlobalDofSection, NULL, "-vertex_global_section_view");CHKERRQ(ierr); 146624121865SAdrian Maldonado PetscFunctionReturn(0); 146724121865SAdrian Maldonado } 14687b6afd5bSHong Zhang 1469f11a936eSBarry Smith /* 1470f11a936eSBarry Smith Add all subnetid for the input vertex v in this process to the btable 1471f11a936eSBarry Smith vertex_subnetid = supportingedge_subnetid 1472f11a936eSBarry Smith */ 1473*9fbee547SJacob Faibussowitsch static inline PetscErrorCode SetSubnetIdLookupBT(DM dm,PetscInt v,PetscInt Nsubnet,PetscBT btable) 1474f11a936eSBarry Smith { 1475f11a936eSBarry Smith PetscErrorCode ierr; 1476f11a936eSBarry Smith PetscInt e,nedges,offset; 1477f11a936eSBarry Smith const PetscInt *edges; 1478f11a936eSBarry Smith DM_Network *newDMnetwork = (DM_Network*)dm->data; 1479f11a936eSBarry Smith DMNetworkComponentHeader header; 1480f11a936eSBarry Smith 1481f11a936eSBarry Smith PetscFunctionBegin; 1482f11a936eSBarry Smith ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 1483f11a936eSBarry Smith ierr = PetscBTMemzero(Nsubnet,btable);CHKERRQ(ierr); 1484f11a936eSBarry Smith for (e=0; e<nedges; e++) { 1485f11a936eSBarry Smith ierr = PetscSectionGetOffset(newDMnetwork->DataSection,edges[e],&offset);CHKERRQ(ierr); 1486f11a936eSBarry Smith header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset); 1487f11a936eSBarry Smith ierr = PetscBTSet(btable,header->subnetid);CHKERRQ(ierr); 1488f11a936eSBarry Smith } 1489f11a936eSBarry Smith PetscFunctionReturn(0); 1490f11a936eSBarry Smith } 1491f11a936eSBarry Smith 14925f2c45f1SShri Abhyankar /*@ 14932bf73ac6SHong Zhang DMNetworkDistribute - Distributes the network and moves associated component data 14945f2c45f1SShri Abhyankar 14955f2c45f1SShri Abhyankar Collective 14965f2c45f1SShri Abhyankar 1497d8d19677SJose E. Roman Input Parameters: 1498d3464fd4SAdrian Maldonado + DM - the DMNetwork object 14992bf73ac6SHong Zhang - overlap - the overlap of partitions, 0 is the default 15005f2c45f1SShri Abhyankar 15015b3f975aSHong Zhang Options Database Key: 15025b3f975aSHong Zhang + -dmnetwork_view - Calls DMView() at the conclusion of DMSetUp() 15035b3f975aSHong Zhang - -dmnetwork_view_distributed - Calls DMView() at the conclusion of DMNetworkDistribute() 15045b3f975aSHong Zhang 15055f2c45f1SShri Abhyankar Notes: 15068b171c8eSHong Zhang Distributes the network with <overlap>-overlapping partitioning of the edges. 15075f2c45f1SShri Abhyankar 15085f2c45f1SShri Abhyankar Level: intermediate 15095f2c45f1SShri Abhyankar 15102bf73ac6SHong Zhang .seealso: DMNetworkCreate() 15115f2c45f1SShri Abhyankar @*/ 1512d3464fd4SAdrian Maldonado PetscErrorCode DMNetworkDistribute(DM *dm,PetscInt overlap) 15135f2c45f1SShri Abhyankar { 1514d3464fd4SAdrian Maldonado MPI_Comm comm; 15155f2c45f1SShri Abhyankar PetscErrorCode ierr; 1516d3464fd4SAdrian Maldonado PetscMPIInt size; 1517d3464fd4SAdrian Maldonado DM_Network *oldDMnetwork = (DM_Network*)((*dm)->data); 1518d3464fd4SAdrian Maldonado DM_Network *newDMnetwork; 1519caf410d2SHong Zhang PetscSF pointsf=NULL; 15205f2c45f1SShri Abhyankar DM newDM; 1521f11a936eSBarry Smith PetscInt j,e,v,offset,*subnetvtx,*subnetedge,Nsubnet,gidx,svtx_idx,nv; 15222bf73ac6SHong Zhang PetscInt to_net,from_net,*svto; 1523f11a936eSBarry Smith PetscBT btable; 152451ac5effSHong Zhang PetscPartitioner part; 1525b9c6e19dSShri Abhyankar DMNetworkComponentHeader header; 15265f2c45f1SShri Abhyankar 15275f2c45f1SShri Abhyankar PetscFunctionBegin; 1528d3464fd4SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)*dm,&comm);CHKERRQ(ierr); 1529ffc4695bSBarry Smith ierr = MPI_Comm_size(comm, &size);CHKERRMPI(ierr); 1530d3464fd4SAdrian Maldonado if (size == 1) PetscFunctionReturn(0); 1531d3464fd4SAdrian Maldonado 15329ace16cdSJacob Faibussowitsch PetscAssertFalse(overlap,PetscObjectComm((PetscObject)dm),PETSC_ERR_SUP,"overlap %D != 0 is not supported yet",overlap); 15335b3f975aSHong Zhang 15342bf73ac6SHong Zhang /* This routine moves the component data to the appropriate processors. It makes use of the DataSection and the componentdataarray to move the component data to appropriate processors and returns a new DataSection and new componentdataarray. */ 1535d3464fd4SAdrian Maldonado ierr = DMNetworkCreate(PetscObjectComm((PetscObject)*dm),&newDM);CHKERRQ(ierr); 15365f2c45f1SShri Abhyankar newDMnetwork = (DM_Network*)newDM->data; 153754dfd506SHong Zhang newDMnetwork->max_comps_registered = oldDMnetwork->max_comps_registered; 153854dfd506SHong Zhang ierr = PetscMalloc1(newDMnetwork->max_comps_registered,&newDMnetwork->component);CHKERRQ(ierr); 153951ac5effSHong Zhang 154051ac5effSHong Zhang /* Enable runtime options for petscpartitioner */ 154151ac5effSHong Zhang ierr = DMPlexGetPartitioner(oldDMnetwork->plex,&part);CHKERRQ(ierr); 154251ac5effSHong Zhang ierr = PetscPartitionerSetFromOptions(part);CHKERRQ(ierr); 154351ac5effSHong Zhang 15442bf73ac6SHong Zhang /* Distribute plex dm */ 154580cf41d5SMatthew G. Knepley ierr = DMPlexDistribute(oldDMnetwork->plex,overlap,&pointsf,&newDMnetwork->plex);CHKERRQ(ierr); 154651ac5effSHong Zhang 15475f2c45f1SShri Abhyankar /* Distribute dof section */ 15482bf73ac6SHong Zhang ierr = PetscSectionCreate(comm,&newDMnetwork->DofSection);CHKERRQ(ierr); 15495f2c45f1SShri Abhyankar ierr = PetscSFDistributeSection(pointsf,oldDMnetwork->DofSection,NULL,newDMnetwork->DofSection);CHKERRQ(ierr); 155051ac5effSHong Zhang 15515f2c45f1SShri Abhyankar /* Distribute data and associated section */ 15522bf73ac6SHong Zhang ierr = PetscSectionCreate(comm,&newDMnetwork->DataSection);CHKERRQ(ierr); 155331da1fc8SHong Zhang ierr = DMPlexDistributeData(newDMnetwork->plex,pointsf,oldDMnetwork->DataSection,MPIU_INT,(void*)oldDMnetwork->componentdataarray,newDMnetwork->DataSection,(void**)&newDMnetwork->componentdataarray);CHKERRQ(ierr); 155424121865SAdrian Maldonado 15555f2c45f1SShri Abhyankar ierr = PetscSectionGetChart(newDMnetwork->DataSection,&newDMnetwork->pStart,&newDMnetwork->pEnd);CHKERRQ(ierr); 15565f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,0, &newDMnetwork->eStart,&newDMnetwork->eEnd);CHKERRQ(ierr); 15575f2c45f1SShri Abhyankar ierr = DMPlexGetHeightStratum(newDMnetwork->plex,1,&newDMnetwork->vStart,&newDMnetwork->vEnd);CHKERRQ(ierr); 15585f2c45f1SShri Abhyankar newDMnetwork->nEdges = newDMnetwork->eEnd - newDMnetwork->eStart; 15596fefedf4SHong Zhang newDMnetwork->nVertices = newDMnetwork->vEnd - newDMnetwork->vStart; 15606fefedf4SHong Zhang newDMnetwork->NVertices = oldDMnetwork->NVertices; 15615f2c45f1SShri Abhyankar newDMnetwork->NEdges = oldDMnetwork->NEdges; 1562f11a936eSBarry Smith newDMnetwork->svtable = oldDMnetwork->svtable; /* global table! */ 15632bf73ac6SHong Zhang oldDMnetwork->svtable = NULL; 156424121865SAdrian Maldonado 15651bb6d2a8SBarry Smith /* Set Dof section as the section for dm */ 156692fd8e1eSJed Brown ierr = DMSetLocalSection(newDMnetwork->plex,newDMnetwork->DofSection);CHKERRQ(ierr); 1567e87a4003SBarry Smith ierr = DMGetGlobalSection(newDMnetwork->plex,&newDMnetwork->GlobalDofSection);CHKERRQ(ierr); 15685f2c45f1SShri Abhyankar 1569b9c6e19dSShri Abhyankar /* Setup subnetwork info in the newDM */ 15702bf73ac6SHong Zhang newDMnetwork->Nsubnet = oldDMnetwork->Nsubnet; 15712bf73ac6SHong Zhang newDMnetwork->Nsvtx = oldDMnetwork->Nsvtx; 15722bf73ac6SHong Zhang oldDMnetwork->Nsvtx = 0; 15732bf73ac6SHong Zhang newDMnetwork->svtx = oldDMnetwork->svtx; /* global vertices! */ 15742bf73ac6SHong Zhang oldDMnetwork->svtx = NULL; 15752bf73ac6SHong Zhang ierr = PetscCalloc1(newDMnetwork->Nsubnet,&newDMnetwork->subnet);CHKERRQ(ierr); 15762bf73ac6SHong Zhang 15772bf73ac6SHong Zhang /* Copy over the global number of vertices and edges in each subnetwork. 15782bf73ac6SHong Zhang Note: these are calculated in DMNetworkLayoutSetUp() 1579b9c6e19dSShri Abhyankar */ 1580f11a936eSBarry Smith Nsubnet = newDMnetwork->Nsubnet; 1581f11a936eSBarry Smith for (j = 0; j < Nsubnet; j++) { 1582b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].Nvtx = oldDMnetwork->subnet[j].Nvtx; 1583b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].Nedge = oldDMnetwork->subnet[j].Nedge; 1584b9c6e19dSShri Abhyankar } 1585b9c6e19dSShri Abhyankar 1586f11a936eSBarry Smith /* Count local nedges for subnetworks */ 1587b9c6e19dSShri Abhyankar for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) { 1588b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 158954dfd506SHong Zhang header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset); 159054dfd506SHong Zhang /* Update pointers */ 159154dfd506SHong Zhang header->size = (PetscInt*)(header + 1); 159254dfd506SHong Zhang header->key = header->size + header->maxcomps; 159354dfd506SHong Zhang header->offset = header->key + header->maxcomps; 159454dfd506SHong Zhang header->nvar = header->offset + header->maxcomps; 159554dfd506SHong Zhang header->offsetvarrel = header->nvar + header->maxcomps; 159654dfd506SHong Zhang 1597b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].nedge++; 1598b9c6e19dSShri Abhyankar } 1599b9c6e19dSShri Abhyankar 1600f11a936eSBarry Smith /* Count local nvtx for subnetworks */ 1601f11a936eSBarry Smith ierr = PetscBTCreate(Nsubnet,&btable);CHKERRQ(ierr); 1602b9c6e19dSShri Abhyankar for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) { 1603b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1604b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 160554dfd506SHong Zhang /* Update pointers */ 160654dfd506SHong Zhang header->size = (PetscInt*)(header + 1); 160754dfd506SHong Zhang header->key = header->size + header->maxcomps; 160854dfd506SHong Zhang header->offset = header->key + header->maxcomps; 160954dfd506SHong Zhang header->nvar = header->offset + header->maxcomps; 161054dfd506SHong Zhang header->offsetvarrel = header->nvar + header->maxcomps; 161154dfd506SHong Zhang 16122bf73ac6SHong Zhang /* shared vertices: use gidx=header->index to check if v is a shared vertex */ 16132bf73ac6SHong Zhang gidx = header->index; 16142bf73ac6SHong Zhang ierr = PetscTableFind(newDMnetwork->svtable,gidx+1,&svtx_idx);CHKERRQ(ierr); 16152bf73ac6SHong Zhang svtx_idx--; 16162bf73ac6SHong Zhang 16172bf73ac6SHong Zhang if (svtx_idx < 0) { /* not a shared vertex */ 1618b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].nvtx++; 16192bf73ac6SHong Zhang } else { /* a shared vertex belongs to more than one subnetworks, it is being counted by multiple subnets */ 1620f11a936eSBarry Smith ierr = SetSubnetIdLookupBT(newDM,v,Nsubnet,btable);CHKERRQ(ierr); 1621f11a936eSBarry Smith 16222bf73ac6SHong Zhang from_net = newDMnetwork->svtx[svtx_idx].sv[0]; 1623f11a936eSBarry Smith if (PetscBTLookup(btable,from_net)) newDMnetwork->subnet[from_net].nvtx++; /* sv is on from_net */ 1624f11a936eSBarry Smith 16252bf73ac6SHong Zhang for (j=1; j<newDMnetwork->svtx[svtx_idx].n; j++) { 16262bf73ac6SHong Zhang svto = newDMnetwork->svtx[svtx_idx].sv + 2*j; 16272bf73ac6SHong Zhang to_net = svto[0]; 1628f11a936eSBarry Smith if (PetscBTLookup(btable,to_net)) newDMnetwork->subnet[to_net].nvtx++; /* sv is on to_net */ 16292bf73ac6SHong Zhang } 16302bf73ac6SHong Zhang } 1631b9c6e19dSShri Abhyankar } 1632b9c6e19dSShri Abhyankar 16332bf73ac6SHong Zhang /* Get total local nvtx for subnetworks */ 16342bf73ac6SHong Zhang nv = 0; 1635f11a936eSBarry Smith for (j=0; j<Nsubnet; j++) nv += newDMnetwork->subnet[j].nvtx; 16362bf73ac6SHong Zhang nv += newDMnetwork->Nsvtx; 1637caf410d2SHong Zhang 16382bf73ac6SHong Zhang /* Now create the vertices and edge arrays for the subnetworks */ 1639f11a936eSBarry Smith ierr = PetscCalloc2(newDMnetwork->nEdges,&subnetedge,nv,&subnetvtx);CHKERRQ(ierr); /* Maps local vertex to local subnetwork's vertex */ 1640f11a936eSBarry Smith newDMnetwork->subnetedge = subnetedge; 16412bf73ac6SHong Zhang newDMnetwork->subnetvtx = subnetvtx; 1642f11a936eSBarry Smith for (j=0; j < newDMnetwork->Nsubnet; j++) { 1643f11a936eSBarry Smith newDMnetwork->subnet[j].edges = subnetedge; 1644f11a936eSBarry Smith subnetedge += newDMnetwork->subnet[j].nedge; 16452bf73ac6SHong Zhang 1646caf410d2SHong Zhang newDMnetwork->subnet[j].vertices = subnetvtx; 1647caf410d2SHong Zhang subnetvtx += newDMnetwork->subnet[j].nvtx; 1648caf410d2SHong Zhang 16492bf73ac6SHong Zhang /* Temporarily setting nvtx and nedge to 0 so we can use them as counters in the below for loop. These get updated when the vertices and edges are added. */ 1650b9c6e19dSShri Abhyankar newDMnetwork->subnet[j].nvtx = newDMnetwork->subnet[j].nedge = 0; 1651b9c6e19dSShri Abhyankar } 16522bf73ac6SHong Zhang newDMnetwork->svertices = subnetvtx; 1653b9c6e19dSShri Abhyankar 16542bf73ac6SHong Zhang /* Set the edges and vertices in each subnetwork */ 1655b9c6e19dSShri Abhyankar for (e = newDMnetwork->eStart; e < newDMnetwork->eEnd; e++) { 1656b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,e,&offset);CHKERRQ(ierr); 1657b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 1658b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].edges[newDMnetwork->subnet[header->subnetid].nedge++] = e; 1659b9c6e19dSShri Abhyankar } 1660b9c6e19dSShri Abhyankar 16612bf73ac6SHong Zhang nv = 0; 1662b9c6e19dSShri Abhyankar for (v = newDMnetwork->vStart; v < newDMnetwork->vEnd; v++) { 1663b9c6e19dSShri Abhyankar ierr = PetscSectionGetOffset(newDMnetwork->DataSection,v,&offset);CHKERRQ(ierr); 1664b9c6e19dSShri Abhyankar header = (DMNetworkComponentHeader)(newDMnetwork->componentdataarray+offset);CHKERRQ(ierr); 16652bf73ac6SHong Zhang 16662bf73ac6SHong Zhang /* coupling vertices: use gidx = header->index to check if v is a coupling vertex */ 16672bf73ac6SHong Zhang ierr = PetscTableFind(newDMnetwork->svtable,header->index+1,&svtx_idx);CHKERRQ(ierr); 16682bf73ac6SHong Zhang svtx_idx--; 16692bf73ac6SHong Zhang if (svtx_idx < 0) { 1670b9c6e19dSShri Abhyankar newDMnetwork->subnet[header->subnetid].vertices[newDMnetwork->subnet[header->subnetid].nvtx++] = v; 16712bf73ac6SHong Zhang } else { /* a shared vertex */ 16722bf73ac6SHong Zhang newDMnetwork->svertices[nv++] = v; 16732bf73ac6SHong Zhang 1674f11a936eSBarry Smith /* add all subnetid for this shared vertex in this process to btable */ 1675f11a936eSBarry Smith ierr = SetSubnetIdLookupBT(newDM,v,Nsubnet,btable);CHKERRQ(ierr); 1676f11a936eSBarry Smith 16772bf73ac6SHong Zhang from_net = newDMnetwork->svtx[svtx_idx].sv[0]; 1678f11a936eSBarry Smith if (PetscBTLookup(btable,from_net)) 16792bf73ac6SHong Zhang newDMnetwork->subnet[from_net].vertices[newDMnetwork->subnet[from_net].nvtx++] = v; 16802bf73ac6SHong Zhang 16812bf73ac6SHong Zhang for (j=1; j<newDMnetwork->svtx[svtx_idx].n; j++) { 16822bf73ac6SHong Zhang svto = newDMnetwork->svtx[svtx_idx].sv + 2*j; 16832bf73ac6SHong Zhang to_net = svto[0]; 1684f11a936eSBarry Smith if (PetscBTLookup(btable,to_net)) 16852bf73ac6SHong Zhang newDMnetwork->subnet[to_net].vertices[newDMnetwork->subnet[to_net].nvtx++] = v; 1686b9c6e19dSShri Abhyankar } 16872bf73ac6SHong Zhang } 16882bf73ac6SHong Zhang } 16892bf73ac6SHong Zhang newDMnetwork->nsvtx = nv; /* num of local shared vertices */ 1690b9c6e19dSShri Abhyankar 1691caf410d2SHong Zhang newDM->setupcalled = (*dm)->setupcalled; 169222bbedd7SHong Zhang newDMnetwork->distributecalled = PETSC_TRUE; 1693caf410d2SHong Zhang 16942bf73ac6SHong Zhang /* Free spaces */ 169524121865SAdrian Maldonado ierr = PetscSFDestroy(&pointsf);CHKERRQ(ierr); 1696d3464fd4SAdrian Maldonado ierr = DMDestroy(dm);CHKERRQ(ierr); 1697f11a936eSBarry Smith ierr = PetscBTDestroy(&btable);CHKERRQ(ierr); 16982bf73ac6SHong Zhang 16995b3f975aSHong Zhang /* View distributed dmnetwork */ 17005b3f975aSHong Zhang ierr = DMViewFromOptions(newDM,NULL,"-dmnetwork_view_distributed");CHKERRQ(ierr); 17015b3f975aSHong Zhang 1702d3464fd4SAdrian Maldonado *dm = newDM; 17035f2c45f1SShri Abhyankar PetscFunctionReturn(0); 17045f2c45f1SShri Abhyankar } 17055f2c45f1SShri Abhyankar 170624121865SAdrian Maldonado /*@C 17072bf73ac6SHong Zhang PetscSFGetSubSF - Returns an SF for a specific subset of points. Leaves are re-numbered to reflect the new ordering 17082bf73ac6SHong Zhang 17092bf73ac6SHong Zhang Collective 171024121865SAdrian Maldonado 171124121865SAdrian Maldonado Input Parameters: 17129dddd249SSatish Balay + mainSF - the original SF structure 171324121865SAdrian Maldonado - map - a ISLocalToGlobal mapping that contains the subset of points 171424121865SAdrian Maldonado 171597bb3fdcSJose E. Roman Output Parameter: 17169dddd249SSatish Balay . subSF - a subset of the mainSF for the desired subset. 1717ee300463SSatish Balay 1718ee300463SSatish Balay Level: intermediate 17197d928bffSSatish Balay @*/ 17209dddd249SSatish Balay PetscErrorCode PetscSFGetSubSF(PetscSF mainsf,ISLocalToGlobalMapping map,PetscSF *subSF) 17212bf73ac6SHong Zhang { 172224121865SAdrian Maldonado PetscErrorCode ierr; 172324121865SAdrian Maldonado PetscInt nroots, nleaves, *ilocal_sub; 172424121865SAdrian Maldonado PetscInt i, *ilocal_map, nroots_sub, nleaves_sub = 0; 172524121865SAdrian Maldonado PetscInt *local_points, *remote_points; 172624121865SAdrian Maldonado PetscSFNode *iremote_sub; 172724121865SAdrian Maldonado const PetscInt *ilocal; 172824121865SAdrian Maldonado const PetscSFNode *iremote; 172924121865SAdrian Maldonado 173024121865SAdrian Maldonado PetscFunctionBegin; 17319dddd249SSatish Balay ierr = PetscSFGetGraph(mainsf,&nroots,&nleaves,&ilocal,&iremote);CHKERRQ(ierr); 173224121865SAdrian Maldonado 173324121865SAdrian Maldonado /* Look for leaves that pertain to the subset of points. Get the local ordering */ 173424121865SAdrian Maldonado ierr = PetscMalloc1(nleaves,&ilocal_map);CHKERRQ(ierr); 173524121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nleaves,ilocal,NULL,ilocal_map);CHKERRQ(ierr); 173624121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 173724121865SAdrian Maldonado if (ilocal_map[i] != -1) nleaves_sub += 1; 173824121865SAdrian Maldonado } 173924121865SAdrian Maldonado /* Re-number ilocal with subset numbering. Need information from roots */ 174024121865SAdrian Maldonado ierr = PetscMalloc2(nroots,&local_points,nroots,&remote_points);CHKERRQ(ierr); 174124121865SAdrian Maldonado for (i = 0; i < nroots; i++) local_points[i] = i; 174224121865SAdrian Maldonado ierr = ISGlobalToLocalMappingApply(map,IS_GTOLM_MASK,nroots,local_points,NULL,local_points);CHKERRQ(ierr); 1743ad227feaSJunchao Zhang ierr = PetscSFBcastBegin(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE);CHKERRQ(ierr); 1744ad227feaSJunchao Zhang ierr = PetscSFBcastEnd(mainsf, MPIU_INT, local_points, remote_points,MPI_REPLACE);CHKERRQ(ierr); 174524121865SAdrian Maldonado /* Fill up graph using local (that is, local to the subset) numbering. */ 17464b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&ilocal_sub);CHKERRQ(ierr); 17474b70a8deSAdrian Maldonado ierr = PetscMalloc1(nleaves_sub,&iremote_sub);CHKERRQ(ierr); 174824121865SAdrian Maldonado nleaves_sub = 0; 174924121865SAdrian Maldonado for (i = 0; i < nleaves; i++) { 175024121865SAdrian Maldonado if (ilocal_map[i] != -1) { 175124121865SAdrian Maldonado ilocal_sub[nleaves_sub] = ilocal_map[i]; 17524b70a8deSAdrian Maldonado iremote_sub[nleaves_sub].rank = iremote[i].rank; 175324121865SAdrian Maldonado iremote_sub[nleaves_sub].index = remote_points[ilocal[i]]; 175424121865SAdrian Maldonado nleaves_sub += 1; 175524121865SAdrian Maldonado } 175624121865SAdrian Maldonado } 175724121865SAdrian Maldonado ierr = PetscFree2(local_points,remote_points);CHKERRQ(ierr); 175824121865SAdrian Maldonado ierr = ISLocalToGlobalMappingGetSize(map,&nroots_sub);CHKERRQ(ierr); 175924121865SAdrian Maldonado 176024121865SAdrian Maldonado /* Create new subSF */ 176124121865SAdrian Maldonado ierr = PetscSFCreate(PETSC_COMM_WORLD,subSF);CHKERRQ(ierr); 176224121865SAdrian Maldonado ierr = PetscSFSetFromOptions(*subSF);CHKERRQ(ierr); 17634b70a8deSAdrian Maldonado ierr = PetscSFSetGraph(*subSF,nroots_sub,nleaves_sub,ilocal_sub,PETSC_OWN_POINTER,iremote_sub,PETSC_COPY_VALUES);CHKERRQ(ierr); 176424121865SAdrian Maldonado ierr = PetscFree(ilocal_map);CHKERRQ(ierr); 17654b70a8deSAdrian Maldonado ierr = PetscFree(iremote_sub);CHKERRQ(ierr); 176624121865SAdrian Maldonado PetscFunctionReturn(0); 176724121865SAdrian Maldonado } 176824121865SAdrian Maldonado 17695f2c45f1SShri Abhyankar /*@C 17705f2c45f1SShri Abhyankar DMNetworkGetSupportingEdges - Return the supporting edges for this vertex point 17715f2c45f1SShri Abhyankar 17725f2c45f1SShri Abhyankar Not Collective 17735f2c45f1SShri Abhyankar 17745f2c45f1SShri Abhyankar Input Parameters: 17752bf73ac6SHong Zhang + dm - the DMNetwork object 17765f2c45f1SShri Abhyankar - p - the vertex point 17775f2c45f1SShri Abhyankar 1778fd292e60Sprj- Output Parameters: 17795f2c45f1SShri Abhyankar + nedges - number of edges connected to this vertex point 17802bf73ac6SHong Zhang - edges - list of edge points 17815f2c45f1SShri Abhyankar 178297bb938eSShri Abhyankar Level: beginner 17835f2c45f1SShri Abhyankar 17845f2c45f1SShri Abhyankar Fortran Notes: 17855f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 17865f2c45f1SShri Abhyankar include petsc.h90 in your code. 17875f2c45f1SShri Abhyankar 17882bf73ac6SHong Zhang .seealso: DMNetworkCreate(), DMNetworkGetConnectedVertices() 17895f2c45f1SShri Abhyankar @*/ 17905f2c45f1SShri Abhyankar PetscErrorCode DMNetworkGetSupportingEdges(DM dm,PetscInt vertex,PetscInt *nedges,const PetscInt *edges[]) 17915f2c45f1SShri Abhyankar { 17925f2c45f1SShri Abhyankar PetscErrorCode ierr; 17935f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 17945f2c45f1SShri Abhyankar 17955f2c45f1SShri Abhyankar PetscFunctionBegin; 17965f2c45f1SShri Abhyankar ierr = DMPlexGetSupportSize(network->plex,vertex,nedges);CHKERRQ(ierr); 1797a9b4a83eSHong Zhang if (edges) {ierr = DMPlexGetSupport(network->plex,vertex,edges);CHKERRQ(ierr);} 17985f2c45f1SShri Abhyankar PetscFunctionReturn(0); 17995f2c45f1SShri Abhyankar } 18005f2c45f1SShri Abhyankar 18015f2c45f1SShri Abhyankar /*@C 1802d842c372SHong Zhang DMNetworkGetConnectedVertices - Return the connected vertices for this edge point 18035f2c45f1SShri Abhyankar 18045f2c45f1SShri Abhyankar Not Collective 18055f2c45f1SShri Abhyankar 18065f2c45f1SShri Abhyankar Input Parameters: 18072bf73ac6SHong Zhang + dm - the DMNetwork object 18085f2c45f1SShri Abhyankar - p - the edge point 18095f2c45f1SShri Abhyankar 1810fd292e60Sprj- Output Parameters: 18115f2c45f1SShri Abhyankar . vertices - vertices connected to this edge 18125f2c45f1SShri Abhyankar 181397bb938eSShri Abhyankar Level: beginner 18145f2c45f1SShri Abhyankar 18155f2c45f1SShri Abhyankar Fortran Notes: 18165f2c45f1SShri Abhyankar Since it returns an array, this routine is only available in Fortran 90, and you must 18175f2c45f1SShri Abhyankar include petsc.h90 in your code. 18185f2c45f1SShri Abhyankar 18192bf73ac6SHong Zhang .seealso: DMNetworkCreate(), DMNetworkGetSupportingEdges() 18205f2c45f1SShri Abhyankar @*/ 1821d842c372SHong Zhang PetscErrorCode DMNetworkGetConnectedVertices(DM dm,PetscInt edge,const PetscInt *vertices[]) 18225f2c45f1SShri Abhyankar { 18235f2c45f1SShri Abhyankar PetscErrorCode ierr; 18245f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 18255f2c45f1SShri Abhyankar 18265f2c45f1SShri Abhyankar PetscFunctionBegin; 18275f2c45f1SShri Abhyankar ierr = DMPlexGetCone(network->plex,edge,vertices);CHKERRQ(ierr); 18285f2c45f1SShri Abhyankar PetscFunctionReturn(0); 18295f2c45f1SShri Abhyankar } 18305f2c45f1SShri Abhyankar 18315f2c45f1SShri Abhyankar /*@ 18322bf73ac6SHong Zhang DMNetworkIsSharedVertex - Returns TRUE if the vertex is shared by subnetworks 18332bf73ac6SHong Zhang 18342bf73ac6SHong Zhang Not Collective 18352bf73ac6SHong Zhang 18362bf73ac6SHong Zhang Input Parameters: 18372bf73ac6SHong Zhang + dm - the DMNetwork object 18382bf73ac6SHong Zhang - p - the vertex point 18392bf73ac6SHong Zhang 18402bf73ac6SHong Zhang Output Parameter: 18412bf73ac6SHong Zhang . flag - TRUE if the vertex is shared by subnetworks 18422bf73ac6SHong Zhang 18432bf73ac6SHong Zhang Level: beginner 18442bf73ac6SHong Zhang 18452bf73ac6SHong Zhang .seealso: DMNetworkAddSharedVertices(), DMNetworkIsGhostVertex() 18462bf73ac6SHong Zhang @*/ 18472bf73ac6SHong Zhang PetscErrorCode DMNetworkIsSharedVertex(DM dm,PetscInt p,PetscBool *flag) 18482bf73ac6SHong Zhang { 18492bf73ac6SHong Zhang PetscErrorCode ierr; 18502bf73ac6SHong Zhang PetscInt i; 18512bf73ac6SHong Zhang 18522bf73ac6SHong Zhang PetscFunctionBegin; 18532bf73ac6SHong Zhang *flag = PETSC_FALSE; 18542bf73ac6SHong Zhang 18552bf73ac6SHong Zhang if (dm->setupcalled) { /* DMNetworkGetGlobalVertexIndex() requires DMSetUp() be called */ 18562bf73ac6SHong Zhang DM_Network *network = (DM_Network*)dm->data; 18572bf73ac6SHong Zhang PetscInt gidx; 18582bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVertexIndex(dm,p,&gidx);CHKERRQ(ierr); 18592bf73ac6SHong Zhang ierr = PetscTableFind(network->svtable,gidx+1,&i);CHKERRQ(ierr); 18602bf73ac6SHong Zhang if (i) *flag = PETSC_TRUE; 18612bf73ac6SHong Zhang } else { /* would be removed? */ 18622bf73ac6SHong Zhang PetscInt nv; 18632bf73ac6SHong Zhang const PetscInt *vtx; 18642bf73ac6SHong Zhang ierr = DMNetworkGetSharedVertices(dm,&nv,&vtx);CHKERRQ(ierr); 18652bf73ac6SHong Zhang for (i=0; i<nv; i++) { 18662bf73ac6SHong Zhang if (p == vtx[i]) { 18672bf73ac6SHong Zhang *flag = PETSC_TRUE; 18682bf73ac6SHong Zhang break; 18692bf73ac6SHong Zhang } 18702bf73ac6SHong Zhang } 18712bf73ac6SHong Zhang } 18722bf73ac6SHong Zhang PetscFunctionReturn(0); 18732bf73ac6SHong Zhang } 18742bf73ac6SHong Zhang 18752bf73ac6SHong Zhang /*@ 18765f2c45f1SShri Abhyankar DMNetworkIsGhostVertex - Returns TRUE if the vertex is a ghost vertex 18775f2c45f1SShri Abhyankar 18785f2c45f1SShri Abhyankar Not Collective 18795f2c45f1SShri Abhyankar 18805f2c45f1SShri Abhyankar Input Parameters: 18812bf73ac6SHong Zhang + dm - the DMNetwork object 1882a2b725a8SWilliam Gropp - p - the vertex point 18835f2c45f1SShri Abhyankar 18845f2c45f1SShri Abhyankar Output Parameter: 18855f2c45f1SShri Abhyankar . isghost - TRUE if the vertex is a ghost point 18865f2c45f1SShri Abhyankar 188797bb938eSShri Abhyankar Level: beginner 18885f2c45f1SShri Abhyankar 18892bf73ac6SHong Zhang .seealso: DMNetworkGetConnectedVertices(), DMNetworkGetVertexRange(), DMNetworkIsSharedVertex() 18905f2c45f1SShri Abhyankar @*/ 18915f2c45f1SShri Abhyankar PetscErrorCode DMNetworkIsGhostVertex(DM dm,PetscInt p,PetscBool *isghost) 18925f2c45f1SShri Abhyankar { 18935f2c45f1SShri Abhyankar PetscErrorCode ierr; 18945f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 18955f2c45f1SShri Abhyankar PetscInt offsetg; 18965f2c45f1SShri Abhyankar PetscSection sectiong; 18975f2c45f1SShri Abhyankar 18985f2c45f1SShri Abhyankar PetscFunctionBegin; 18995f2c45f1SShri Abhyankar *isghost = PETSC_FALSE; 1900e87a4003SBarry Smith ierr = DMGetGlobalSection(network->plex,§iong);CHKERRQ(ierr); 19015f2c45f1SShri Abhyankar ierr = PetscSectionGetOffset(sectiong,p,&offsetg);CHKERRQ(ierr); 19025f2c45f1SShri Abhyankar if (offsetg < 0) *isghost = PETSC_TRUE; 19035f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19045f2c45f1SShri Abhyankar } 19055f2c45f1SShri Abhyankar 19065f2c45f1SShri Abhyankar PetscErrorCode DMSetUp_Network(DM dm) 19075f2c45f1SShri Abhyankar { 19085f2c45f1SShri Abhyankar PetscErrorCode ierr; 19095f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 19105f2c45f1SShri Abhyankar 19115f2c45f1SShri Abhyankar PetscFunctionBegin; 19125f2c45f1SShri Abhyankar ierr = DMNetworkComponentSetUp(dm);CHKERRQ(ierr); 19135f2c45f1SShri Abhyankar ierr = DMNetworkVariablesSetUp(dm);CHKERRQ(ierr); 19145f2c45f1SShri Abhyankar 191592fd8e1eSJed Brown ierr = DMSetLocalSection(network->plex,network->DofSection);CHKERRQ(ierr); 1916e87a4003SBarry Smith ierr = DMGetGlobalSection(network->plex,&network->GlobalDofSection);CHKERRQ(ierr); 19179e1d080bSHong Zhang 19189e1d080bSHong Zhang dm->setupcalled = PETSC_TRUE; 19195b3f975aSHong Zhang 19205b3f975aSHong Zhang /* View dmnetwork */ 19215b3f975aSHong Zhang ierr = DMViewFromOptions(dm,NULL,"-dmnetwork_view");CHKERRQ(ierr); 19225f2c45f1SShri Abhyankar PetscFunctionReturn(0); 19235f2c45f1SShri Abhyankar } 19245f2c45f1SShri Abhyankar 19251ad426b7SHong Zhang /*@ 192617df6e9eSHong Zhang DMNetworkHasJacobian - Sets global flag for using user's sub Jacobian matrices 19271ad426b7SHong Zhang -- replaced by DMNetworkSetOption(network,userjacobian,PETSC_TURE)? 19281ad426b7SHong Zhang 19291ad426b7SHong Zhang Collective 19301ad426b7SHong Zhang 19311ad426b7SHong Zhang Input Parameters: 19322bf73ac6SHong Zhang + dm - the DMNetwork object 193383b2e829SHong Zhang . eflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for edges 193483b2e829SHong Zhang - vflg - turn the option on (PETSC_TRUE) or off (PETSC_FALSE) if user provides Jacobian for vertices 19351ad426b7SHong Zhang 19361ad426b7SHong Zhang Level: intermediate 19371ad426b7SHong Zhang 19381ad426b7SHong Zhang @*/ 193983b2e829SHong Zhang PetscErrorCode DMNetworkHasJacobian(DM dm,PetscBool eflg,PetscBool vflg) 19401ad426b7SHong Zhang { 19411ad426b7SHong Zhang DM_Network *network=(DM_Network*)dm->data; 19428675203cSHong Zhang PetscErrorCode ierr; 194366b4e0ffSHong Zhang PetscInt nVertices = network->nVertices; 19441ad426b7SHong Zhang 19451ad426b7SHong Zhang PetscFunctionBegin; 194683b2e829SHong Zhang network->userEdgeJacobian = eflg; 194783b2e829SHong Zhang network->userVertexJacobian = vflg; 19488675203cSHong Zhang 19498675203cSHong Zhang if (eflg && !network->Je) { 19508675203cSHong Zhang ierr = PetscCalloc1(3*network->nEdges,&network->Je);CHKERRQ(ierr); 19518675203cSHong Zhang } 19528675203cSHong Zhang 195366b4e0ffSHong Zhang if (vflg && !network->Jv && nVertices) { 19548675203cSHong Zhang PetscInt i,*vptr,nedges,vStart=network->vStart; 195566b4e0ffSHong Zhang PetscInt nedges_total; 19568675203cSHong Zhang const PetscInt *edges; 19578675203cSHong Zhang 19588675203cSHong Zhang /* count nvertex_total */ 19598675203cSHong Zhang nedges_total = 0; 19608675203cSHong Zhang ierr = PetscMalloc1(nVertices+1,&vptr);CHKERRQ(ierr); 19618675203cSHong Zhang 19628675203cSHong Zhang vptr[0] = 0; 19638675203cSHong Zhang for (i=0; i<nVertices; i++) { 19648675203cSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,i+vStart,&nedges,&edges);CHKERRQ(ierr); 19658675203cSHong Zhang nedges_total += nedges; 19668675203cSHong Zhang vptr[i+1] = vptr[i] + 2*nedges + 1; 19678675203cSHong Zhang } 19688675203cSHong Zhang 19698675203cSHong Zhang ierr = PetscCalloc1(2*nedges_total+nVertices,&network->Jv);CHKERRQ(ierr); 19708675203cSHong Zhang network->Jvptr = vptr; 19718675203cSHong Zhang } 19721ad426b7SHong Zhang PetscFunctionReturn(0); 19731ad426b7SHong Zhang } 19741ad426b7SHong Zhang 19751ad426b7SHong Zhang /*@ 197683b2e829SHong Zhang DMNetworkEdgeSetMatrix - Sets user-provided Jacobian matrices for this edge to the network 197783b2e829SHong Zhang 197883b2e829SHong Zhang Not Collective 197983b2e829SHong Zhang 198083b2e829SHong Zhang Input Parameters: 19812bf73ac6SHong Zhang + dm - the DMNetwork object 198283b2e829SHong Zhang . p - the edge point 19833e97b6e8SHong Zhang - J - array (size = 3) of Jacobian submatrices for this edge point: 19843e97b6e8SHong Zhang J[0]: this edge 1985d842c372SHong Zhang J[1] and J[2]: connected vertices, obtained by calling DMNetworkGetConnectedVertices() 198683b2e829SHong Zhang 198797bb938eSShri Abhyankar Level: advanced 198883b2e829SHong Zhang 19892bf73ac6SHong Zhang .seealso: DMNetworkVertexSetMatrix() 199083b2e829SHong Zhang @*/ 199183b2e829SHong Zhang PetscErrorCode DMNetworkEdgeSetMatrix(DM dm,PetscInt p,Mat J[]) 199283b2e829SHong Zhang { 199383b2e829SHong Zhang DM_Network *network=(DM_Network*)dm->data; 199483b2e829SHong Zhang 199583b2e829SHong Zhang PetscFunctionBegin; 19969ace16cdSJacob Faibussowitsch PetscAssertFalse(!network->Je,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkEdgeSetMatrix"); 19978675203cSHong Zhang 19988675203cSHong Zhang if (J) { 1999883e35e8SHong Zhang network->Je[3*p] = J[0]; 2000883e35e8SHong Zhang network->Je[3*p+1] = J[1]; 2001883e35e8SHong Zhang network->Je[3*p+2] = J[2]; 20028675203cSHong Zhang } 200383b2e829SHong Zhang PetscFunctionReturn(0); 200483b2e829SHong Zhang } 200583b2e829SHong Zhang 200683b2e829SHong Zhang /*@ 200776ddfea5SHong Zhang DMNetworkVertexSetMatrix - Sets user-provided Jacobian matrix for this vertex to the network 20081ad426b7SHong Zhang 20091ad426b7SHong Zhang Not Collective 20101ad426b7SHong Zhang 20111ad426b7SHong Zhang Input Parameters: 20121ad426b7SHong Zhang + dm - The DMNetwork object 20131ad426b7SHong Zhang . p - the vertex point 20143e97b6e8SHong Zhang - J - array of Jacobian (size = 2*(num of supporting edges) + 1) submatrices for this vertex point: 20153e97b6e8SHong Zhang J[0]: this vertex 20163e97b6e8SHong Zhang J[1+2*i]: i-th supporting edge 20173e97b6e8SHong Zhang J[1+2*i+1]: i-th connected vertex 20181ad426b7SHong Zhang 201997bb938eSShri Abhyankar Level: advanced 20201ad426b7SHong Zhang 20212bf73ac6SHong Zhang .seealso: DMNetworkEdgeSetMatrix() 20221ad426b7SHong Zhang @*/ 2023883e35e8SHong Zhang PetscErrorCode DMNetworkVertexSetMatrix(DM dm,PetscInt p,Mat J[]) 20245f2c45f1SShri Abhyankar { 20255f2c45f1SShri Abhyankar PetscErrorCode ierr; 20265f2c45f1SShri Abhyankar DM_Network *network=(DM_Network*)dm->data; 20278675203cSHong Zhang PetscInt i,*vptr,nedges,vStart=network->vStart; 2028883e35e8SHong Zhang const PetscInt *edges; 20295f2c45f1SShri Abhyankar 20305f2c45f1SShri Abhyankar PetscFunctionBegin; 20319ace16cdSJacob Faibussowitsch PetscAssertFalse(!network->Jv,PetscObjectComm((PetscObject)dm),PETSC_ERR_ORDER,"Must call DMNetworkHasJacobian() collectively before calling DMNetworkVertexSetMatrix"); 2032883e35e8SHong Zhang 20338675203cSHong Zhang if (J) { 2034883e35e8SHong Zhang vptr = network->Jvptr; 20353e97b6e8SHong Zhang network->Jv[vptr[p-vStart]] = J[0]; /* Set Jacobian for this vertex */ 20363e97b6e8SHong Zhang 20373e97b6e8SHong Zhang /* Set Jacobian for each supporting edge and connected vertex */ 2038883e35e8SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,p,&nedges,&edges);CHKERRQ(ierr); 2039883e35e8SHong Zhang for (i=1; i<=2*nedges; i++) network->Jv[vptr[p-vStart]+i] = J[i]; 20408675203cSHong Zhang } 2041883e35e8SHong Zhang PetscFunctionReturn(0); 2042883e35e8SHong Zhang } 2043883e35e8SHong Zhang 2044*9fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 20455cf7da58SHong Zhang { 20465cf7da58SHong Zhang PetscErrorCode ierr; 20475cf7da58SHong Zhang PetscInt j; 20485cf7da58SHong Zhang PetscScalar val=(PetscScalar)ncols; 20495cf7da58SHong Zhang 20505cf7da58SHong Zhang PetscFunctionBegin; 20515cf7da58SHong Zhang if (!ghost) { 20525cf7da58SHong Zhang for (j=0; j<nrows; j++) { 20535cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 20545cf7da58SHong Zhang } 20555cf7da58SHong Zhang } else { 20565cf7da58SHong Zhang for (j=0; j<nrows; j++) { 20575cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 20585cf7da58SHong Zhang } 20595cf7da58SHong Zhang } 20605cf7da58SHong Zhang PetscFunctionReturn(0); 20615cf7da58SHong Zhang } 20625cf7da58SHong Zhang 2063*9fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 20645cf7da58SHong Zhang { 20655cf7da58SHong Zhang PetscErrorCode ierr; 20665cf7da58SHong Zhang PetscInt j,ncols_u; 20675cf7da58SHong Zhang PetscScalar val; 20685cf7da58SHong Zhang 20695cf7da58SHong Zhang PetscFunctionBegin; 20705cf7da58SHong Zhang if (!ghost) { 20715cf7da58SHong Zhang for (j=0; j<nrows; j++) { 20725cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 20735cf7da58SHong Zhang val = (PetscScalar)ncols_u; 20745cf7da58SHong Zhang ierr = VecSetValues(vdnz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 20755cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 20765cf7da58SHong Zhang } 20775cf7da58SHong Zhang } else { 20785cf7da58SHong Zhang for (j=0; j<nrows; j++) { 20795cf7da58SHong Zhang ierr = MatGetRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 20805cf7da58SHong Zhang val = (PetscScalar)ncols_u; 20815cf7da58SHong Zhang ierr = VecSetValues(vonz,1,&rows[j],&val,ADD_VALUES);CHKERRQ(ierr); 20825cf7da58SHong Zhang ierr = MatRestoreRow(Ju,j,&ncols_u,NULL,NULL);CHKERRQ(ierr); 20835cf7da58SHong Zhang } 20845cf7da58SHong Zhang } 20855cf7da58SHong Zhang PetscFunctionReturn(0); 20865cf7da58SHong Zhang } 20875cf7da58SHong Zhang 2088*9fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetPreallocationblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscBool ghost,Vec vdnz,Vec vonz) 20895cf7da58SHong Zhang { 20905cf7da58SHong Zhang PetscErrorCode ierr; 20915cf7da58SHong Zhang 20925cf7da58SHong Zhang PetscFunctionBegin; 20935cf7da58SHong Zhang if (Ju) { 20945cf7da58SHong Zhang ierr = MatSetPreallocationUserblock_private(Ju,nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 20955cf7da58SHong Zhang } else { 20965cf7da58SHong Zhang ierr = MatSetPreallocationDenseblock_private(nrows,rows,ncols,ghost,vdnz,vonz);CHKERRQ(ierr); 20975cf7da58SHong Zhang } 20985cf7da58SHong Zhang PetscFunctionReturn(0); 20995cf7da58SHong Zhang } 21005cf7da58SHong Zhang 2101*9fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetDenseblock_private(PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 2102883e35e8SHong Zhang { 2103883e35e8SHong Zhang PetscErrorCode ierr; 2104883e35e8SHong Zhang PetscInt j,*cols; 2105883e35e8SHong Zhang PetscScalar *zeros; 2106883e35e8SHong Zhang 2107883e35e8SHong Zhang PetscFunctionBegin; 2108883e35e8SHong Zhang ierr = PetscCalloc2(ncols,&cols,nrows*ncols,&zeros);CHKERRQ(ierr); 2109883e35e8SHong Zhang for (j=0; j<ncols; j++) cols[j] = j+ cstart; 2110883e35e8SHong Zhang ierr = MatSetValues(*J,nrows,rows,ncols,cols,zeros,INSERT_VALUES);CHKERRQ(ierr); 2111883e35e8SHong Zhang ierr = PetscFree2(cols,zeros);CHKERRQ(ierr); 21121ad426b7SHong Zhang PetscFunctionReturn(0); 21131ad426b7SHong Zhang } 2114a4e85ca8SHong Zhang 2115*9fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetUserblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 21163e97b6e8SHong Zhang { 21173e97b6e8SHong Zhang PetscErrorCode ierr; 21183e97b6e8SHong Zhang PetscInt j,M,N,row,col,ncols_u; 21193e97b6e8SHong Zhang const PetscInt *cols; 21203e97b6e8SHong Zhang PetscScalar zero=0.0; 21213e97b6e8SHong Zhang 21223e97b6e8SHong Zhang PetscFunctionBegin; 21233e97b6e8SHong Zhang ierr = MatGetSize(Ju,&M,&N);CHKERRQ(ierr); 21249ace16cdSJacob Faibussowitsch PetscAssertFalse(nrows != M || ncols != N,PetscObjectComm((PetscObject)Ju),PETSC_ERR_USER,"%D by %D must equal %D by %D",nrows,ncols,M,N); 21253e97b6e8SHong Zhang 21263e97b6e8SHong Zhang for (row=0; row<nrows; row++) { 21273e97b6e8SHong Zhang ierr = MatGetRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 21283e97b6e8SHong Zhang for (j=0; j<ncols_u; j++) { 21293e97b6e8SHong Zhang col = cols[j] + cstart; 21303e97b6e8SHong Zhang ierr = MatSetValues(*J,1,&rows[row],1,&col,&zero,INSERT_VALUES);CHKERRQ(ierr); 21313e97b6e8SHong Zhang } 21323e97b6e8SHong Zhang ierr = MatRestoreRow(Ju,row,&ncols_u,&cols,NULL);CHKERRQ(ierr); 21333e97b6e8SHong Zhang } 21343e97b6e8SHong Zhang PetscFunctionReturn(0); 21353e97b6e8SHong Zhang } 21361ad426b7SHong Zhang 2137*9fbee547SJacob Faibussowitsch static inline PetscErrorCode MatSetblock_private(Mat Ju,PetscInt nrows,PetscInt *rows,PetscInt ncols,PetscInt cstart,Mat *J) 2138a4e85ca8SHong Zhang { 2139a4e85ca8SHong Zhang PetscErrorCode ierr; 2140f4431b8cSHong Zhang 2141a4e85ca8SHong Zhang PetscFunctionBegin; 2142a4e85ca8SHong Zhang if (Ju) { 2143a4e85ca8SHong Zhang ierr = MatSetUserblock_private(Ju,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 2144a4e85ca8SHong Zhang } else { 2145a4e85ca8SHong Zhang ierr = MatSetDenseblock_private(nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 2146a4e85ca8SHong Zhang } 2147a4e85ca8SHong Zhang PetscFunctionReturn(0); 2148a4e85ca8SHong Zhang } 2149a4e85ca8SHong Zhang 215024121865SAdrian 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. 215124121865SAdrian Maldonado */ 215224121865SAdrian Maldonado PetscErrorCode CreateSubGlobalToLocalMapping_private(PetscSection globalsec, PetscSection localsec, ISLocalToGlobalMapping *ltog) 215324121865SAdrian Maldonado { 215424121865SAdrian Maldonado PetscErrorCode ierr; 215524121865SAdrian Maldonado PetscInt i,size,dof; 215624121865SAdrian Maldonado PetscInt *glob2loc; 215724121865SAdrian Maldonado 215824121865SAdrian Maldonado PetscFunctionBegin; 215924121865SAdrian Maldonado ierr = PetscSectionGetStorageSize(localsec,&size);CHKERRQ(ierr); 216024121865SAdrian Maldonado ierr = PetscMalloc1(size,&glob2loc);CHKERRQ(ierr); 216124121865SAdrian Maldonado 216224121865SAdrian Maldonado for (i = 0; i < size; i++) { 216324121865SAdrian Maldonado ierr = PetscSectionGetOffset(globalsec,i,&dof);CHKERRQ(ierr); 216424121865SAdrian Maldonado dof = (dof >= 0) ? dof : -(dof + 1); 216524121865SAdrian Maldonado glob2loc[i] = dof; 216624121865SAdrian Maldonado } 216724121865SAdrian Maldonado 216824121865SAdrian Maldonado ierr = ISLocalToGlobalMappingCreate(PETSC_COMM_WORLD,1,size,glob2loc,PETSC_OWN_POINTER,ltog);CHKERRQ(ierr); 216924121865SAdrian Maldonado #if 0 217024121865SAdrian Maldonado ierr = PetscIntView(size,glob2loc,PETSC_VIEWER_STDOUT_WORLD);CHKERRQ(ierr); 217124121865SAdrian Maldonado #endif 217224121865SAdrian Maldonado PetscFunctionReturn(0); 217324121865SAdrian Maldonado } 217424121865SAdrian Maldonado 217501ad2aeeSHong Zhang #include <petsc/private/matimpl.h> 21769e1d080bSHong Zhang 21779e1d080bSHong Zhang PetscErrorCode DMCreateMatrix_Network_Nest(DM dm,Mat *J) 21781ad426b7SHong Zhang { 21791ad426b7SHong Zhang PetscErrorCode ierr; 21801ad426b7SHong Zhang DM_Network *network = (DM_Network*)dm->data; 218124121865SAdrian Maldonado PetscInt eDof,vDof; 218224121865SAdrian Maldonado Mat j11,j12,j21,j22,bA[2][2]; 21839e1d080bSHong Zhang MPI_Comm comm; 218424121865SAdrian Maldonado ISLocalToGlobalMapping eISMap,vISMap; 218524121865SAdrian Maldonado 21869e1d080bSHong Zhang PetscFunctionBegin; 218724121865SAdrian Maldonado ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 218824121865SAdrian Maldonado 218924121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->edge.GlobalDofSection,&eDof);CHKERRQ(ierr); 219024121865SAdrian Maldonado ierr = PetscSectionGetConstrainedStorageSize(network->vertex.GlobalDofSection,&vDof);CHKERRQ(ierr); 219124121865SAdrian Maldonado 219201ad2aeeSHong Zhang ierr = MatCreate(comm, &j11);CHKERRQ(ierr); 219324121865SAdrian Maldonado ierr = MatSetSizes(j11, eDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 219424121865SAdrian Maldonado ierr = MatSetType(j11, MATMPIAIJ);CHKERRQ(ierr); 219524121865SAdrian Maldonado 219601ad2aeeSHong Zhang ierr = MatCreate(comm, &j12);CHKERRQ(ierr); 219724121865SAdrian Maldonado ierr = MatSetSizes(j12, eDof, vDof, PETSC_DETERMINE ,PETSC_DETERMINE);CHKERRQ(ierr); 219824121865SAdrian Maldonado ierr = MatSetType(j12, MATMPIAIJ);CHKERRQ(ierr); 219924121865SAdrian Maldonado 220001ad2aeeSHong Zhang ierr = MatCreate(comm, &j21);CHKERRQ(ierr); 220124121865SAdrian Maldonado ierr = MatSetSizes(j21, vDof, eDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 220224121865SAdrian Maldonado ierr = MatSetType(j21, MATMPIAIJ);CHKERRQ(ierr); 220324121865SAdrian Maldonado 220401ad2aeeSHong Zhang ierr = MatCreate(comm, &j22);CHKERRQ(ierr); 220524121865SAdrian Maldonado ierr = MatSetSizes(j22, vDof, vDof, PETSC_DETERMINE, PETSC_DETERMINE);CHKERRQ(ierr); 220624121865SAdrian Maldonado ierr = MatSetType(j22, MATMPIAIJ);CHKERRQ(ierr); 220724121865SAdrian Maldonado 22083f6a6bdaSHong Zhang bA[0][0] = j11; 22093f6a6bdaSHong Zhang bA[0][1] = j12; 22103f6a6bdaSHong Zhang bA[1][0] = j21; 22113f6a6bdaSHong Zhang bA[1][1] = j22; 221224121865SAdrian Maldonado 221324121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->edge.GlobalDofSection,network->edge.DofSection,&eISMap);CHKERRQ(ierr); 221424121865SAdrian Maldonado ierr = CreateSubGlobalToLocalMapping_private(network->vertex.GlobalDofSection,network->vertex.DofSection,&vISMap);CHKERRQ(ierr); 221524121865SAdrian Maldonado 221624121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j11,eISMap,eISMap);CHKERRQ(ierr); 221724121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j12,eISMap,vISMap);CHKERRQ(ierr); 221824121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j21,vISMap,eISMap);CHKERRQ(ierr); 221924121865SAdrian Maldonado ierr = MatSetLocalToGlobalMapping(j22,vISMap,vISMap);CHKERRQ(ierr); 222024121865SAdrian Maldonado 222124121865SAdrian Maldonado ierr = MatSetUp(j11);CHKERRQ(ierr); 222224121865SAdrian Maldonado ierr = MatSetUp(j12);CHKERRQ(ierr); 222324121865SAdrian Maldonado ierr = MatSetUp(j21);CHKERRQ(ierr); 222424121865SAdrian Maldonado ierr = MatSetUp(j22);CHKERRQ(ierr); 222524121865SAdrian Maldonado 222601ad2aeeSHong Zhang ierr = MatCreateNest(comm,2,NULL,2,NULL,&bA[0][0],J);CHKERRQ(ierr); 222724121865SAdrian Maldonado ierr = MatSetUp(*J);CHKERRQ(ierr); 222824121865SAdrian Maldonado ierr = MatNestSetVecType(*J,VECNEST);CHKERRQ(ierr); 222924121865SAdrian Maldonado ierr = MatDestroy(&j11);CHKERRQ(ierr); 223024121865SAdrian Maldonado ierr = MatDestroy(&j12);CHKERRQ(ierr); 223124121865SAdrian Maldonado ierr = MatDestroy(&j21);CHKERRQ(ierr); 223224121865SAdrian Maldonado ierr = MatDestroy(&j22);CHKERRQ(ierr); 223324121865SAdrian Maldonado 223424121865SAdrian Maldonado ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 223524121865SAdrian Maldonado ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 22369e1d080bSHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 223724121865SAdrian Maldonado 223824121865SAdrian Maldonado /* Free structures */ 223924121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&eISMap);CHKERRQ(ierr); 224024121865SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&vISMap);CHKERRQ(ierr); 224124121865SAdrian Maldonado PetscFunctionReturn(0); 22429e1d080bSHong Zhang } 22439e1d080bSHong Zhang 22449e1d080bSHong Zhang PetscErrorCode DMCreateMatrix_Network(DM dm,Mat *J) 22459e1d080bSHong Zhang { 22469e1d080bSHong Zhang PetscErrorCode ierr; 22479e1d080bSHong Zhang DM_Network *network = (DM_Network*)dm->data; 22489e1d080bSHong Zhang PetscInt eStart,eEnd,vStart,vEnd,rstart,nrows,*rows,localSize; 22499e1d080bSHong Zhang PetscInt cstart,ncols,j,e,v; 22509e1d080bSHong Zhang PetscBool ghost,ghost_vc,ghost2,isNest; 22519e1d080bSHong Zhang Mat Juser; 22529e1d080bSHong Zhang PetscSection sectionGlobal; 22539e1d080bSHong Zhang PetscInt nedges,*vptr=NULL,vc,*rows_v; /* suppress maybe-uninitialized warning */ 22549e1d080bSHong Zhang const PetscInt *edges,*cone; 22559e1d080bSHong Zhang MPI_Comm comm; 22569e1d080bSHong Zhang MatType mtype; 22579e1d080bSHong Zhang Vec vd_nz,vo_nz; 22589e1d080bSHong Zhang PetscInt *dnnz,*onnz; 22599e1d080bSHong Zhang PetscScalar *vdnz,*vonz; 22609e1d080bSHong Zhang 22619e1d080bSHong Zhang PetscFunctionBegin; 22629e1d080bSHong Zhang mtype = dm->mattype; 22639e1d080bSHong Zhang ierr = PetscStrcmp(mtype,MATNEST,&isNest);CHKERRQ(ierr); 22649e1d080bSHong Zhang if (isNest) { 22659e1d080bSHong Zhang ierr = DMCreateMatrix_Network_Nest(dm,J);CHKERRQ(ierr); 2266c6b011d8SStefano Zampini ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 22679e1d080bSHong Zhang PetscFunctionReturn(0); 22689e1d080bSHong Zhang } 22699e1d080bSHong Zhang 22709e1d080bSHong Zhang if (!network->userEdgeJacobian && !network->userVertexJacobian) { 2271a4e85ca8SHong Zhang /* user does not provide Jacobian blocks */ 22729e1d080bSHong Zhang ierr = DMCreateMatrix_Plex(network->plex,J);CHKERRQ(ierr); 2273bfbc38dcSHong Zhang ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 22741ad426b7SHong Zhang PetscFunctionReturn(0); 22751ad426b7SHong Zhang } 22761ad426b7SHong Zhang 2277bfbc38dcSHong Zhang ierr = MatCreate(PetscObjectComm((PetscObject)dm),J);CHKERRQ(ierr); 2278e87a4003SBarry Smith ierr = DMGetGlobalSection(network->plex,§ionGlobal);CHKERRQ(ierr); 2279bfbc38dcSHong Zhang ierr = PetscSectionGetConstrainedStorageSize(sectionGlobal,&localSize);CHKERRQ(ierr); 2280bfbc38dcSHong Zhang ierr = MatSetSizes(*J,localSize,localSize,PETSC_DETERMINE,PETSC_DETERMINE);CHKERRQ(ierr); 22812a945128SHong Zhang 22822a945128SHong Zhang ierr = MatSetType(*J,MATAIJ);CHKERRQ(ierr); 22832a945128SHong Zhang ierr = MatSetFromOptions(*J);CHKERRQ(ierr); 228489898e50SHong Zhang 228589898e50SHong Zhang /* (1) Set matrix preallocation */ 228689898e50SHong Zhang /*------------------------------*/ 2287840c2264SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2288840c2264SHong Zhang ierr = VecCreate(comm,&vd_nz);CHKERRQ(ierr); 2289840c2264SHong Zhang ierr = VecSetSizes(vd_nz,localSize,PETSC_DECIDE);CHKERRQ(ierr); 2290840c2264SHong Zhang ierr = VecSetFromOptions(vd_nz);CHKERRQ(ierr); 2291840c2264SHong Zhang ierr = VecSet(vd_nz,0.0);CHKERRQ(ierr); 2292840c2264SHong Zhang ierr = VecDuplicate(vd_nz,&vo_nz);CHKERRQ(ierr); 2293840c2264SHong Zhang 229489898e50SHong Zhang /* Set preallocation for edges */ 229589898e50SHong Zhang /*-----------------------------*/ 2296840c2264SHong Zhang ierr = DMNetworkGetEdgeRange(dm,&eStart,&eEnd);CHKERRQ(ierr); 2297840c2264SHong Zhang 2298bdcb62a2SHong Zhang ierr = PetscMalloc1(localSize,&rows);CHKERRQ(ierr); 2299840c2264SHong Zhang for (e=eStart; e<eEnd; e++) { 2300840c2264SHong Zhang /* Get row indices */ 23012bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart);CHKERRQ(ierr); 23022bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,e,&nrows);CHKERRQ(ierr); 2303840c2264SHong Zhang if (nrows) { 2304840c2264SHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 2305840c2264SHong Zhang 2306a5b23f4aSJose E. Roman /* Set preallocation for connected vertices */ 2307d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 2308840c2264SHong Zhang for (v=0; v<2; v++) { 23092bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,cone[v],&ncols);CHKERRQ(ierr); 2310840c2264SHong Zhang 23118675203cSHong Zhang if (network->Je) { 2312840c2264SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 23138675203cSHong Zhang } else Juser = NULL; 2314840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,cone[v],&ghost);CHKERRQ(ierr); 23155cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 2316840c2264SHong Zhang } 2317840c2264SHong Zhang 231889898e50SHong Zhang /* Set preallocation for edge self */ 2319840c2264SHong Zhang cstart = rstart; 23208675203cSHong Zhang if (network->Je) { 2321840c2264SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 23228675203cSHong Zhang } else Juser = NULL; 23235cf7da58SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 2324840c2264SHong Zhang } 2325840c2264SHong Zhang } 2326840c2264SHong Zhang 232789898e50SHong Zhang /* Set preallocation for vertices */ 232889898e50SHong Zhang /*--------------------------------*/ 2329840c2264SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vStart,&vEnd);CHKERRQ(ierr); 23308675203cSHong Zhang if (vEnd - vStart) vptr = network->Jvptr; 2331840c2264SHong Zhang 2332840c2264SHong Zhang for (v=vStart; v<vEnd; v++) { 2333840c2264SHong Zhang /* Get row indices */ 23342bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart);CHKERRQ(ierr); 23352bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,v,&nrows);CHKERRQ(ierr); 2336840c2264SHong Zhang if (!nrows) continue; 2337840c2264SHong Zhang 2338bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 2339bdcb62a2SHong Zhang if (ghost) { 2340bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 2341bdcb62a2SHong Zhang } else { 2342bdcb62a2SHong Zhang rows_v = rows; 2343bdcb62a2SHong Zhang } 2344bdcb62a2SHong Zhang 2345bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 2346840c2264SHong Zhang 2347840c2264SHong Zhang /* Get supporting edges and connected vertices */ 2348840c2264SHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 2349840c2264SHong Zhang 2350840c2264SHong Zhang for (e=0; e<nedges; e++) { 2351840c2264SHong Zhang /* Supporting edges */ 23522bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart);CHKERRQ(ierr); 23532bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,edges[e],&ncols);CHKERRQ(ierr); 2354840c2264SHong Zhang 23558675203cSHong Zhang if (network->Jv) { 2356840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 23578675203cSHong Zhang } else Juser = NULL; 2358bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost,vd_nz,vo_nz);CHKERRQ(ierr); 2359840c2264SHong Zhang 2360840c2264SHong Zhang /* Connected vertices */ 2361d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 2362840c2264SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 2363840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,vc,&ghost_vc);CHKERRQ(ierr); 2364840c2264SHong Zhang 23652bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,vc,&ncols);CHKERRQ(ierr); 2366840c2264SHong Zhang 23678675203cSHong Zhang if (network->Jv) { 2368840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 23698675203cSHong Zhang } else Juser = NULL; 2370e102a522SHong Zhang if (ghost_vc||ghost) { 2371e102a522SHong Zhang ghost2 = PETSC_TRUE; 2372e102a522SHong Zhang } else { 2373e102a522SHong Zhang ghost2 = PETSC_FALSE; 2374e102a522SHong Zhang } 2375e102a522SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,ncols,ghost2,vd_nz,vo_nz);CHKERRQ(ierr); 2376840c2264SHong Zhang } 2377840c2264SHong Zhang 237889898e50SHong Zhang /* Set preallocation for vertex self */ 2379840c2264SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 2380840c2264SHong Zhang if (!ghost) { 23812bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart);CHKERRQ(ierr); 23828675203cSHong Zhang if (network->Jv) { 2383840c2264SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 23848675203cSHong Zhang } else Juser = NULL; 2385bdcb62a2SHong Zhang ierr = MatSetPreallocationblock_private(Juser,nrows,rows_v,nrows,PETSC_FALSE,vd_nz,vo_nz);CHKERRQ(ierr); 2386840c2264SHong Zhang } 2387bdcb62a2SHong Zhang if (ghost) { 2388bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 2389bdcb62a2SHong Zhang } 2390840c2264SHong Zhang } 2391840c2264SHong Zhang 2392840c2264SHong Zhang ierr = VecAssemblyBegin(vd_nz);CHKERRQ(ierr); 2393840c2264SHong Zhang ierr = VecAssemblyBegin(vo_nz);CHKERRQ(ierr); 23945cf7da58SHong Zhang 23955cf7da58SHong Zhang ierr = PetscMalloc2(localSize,&dnnz,localSize,&onnz);CHKERRQ(ierr); 23965cf7da58SHong Zhang 23975cf7da58SHong Zhang ierr = VecAssemblyEnd(vd_nz);CHKERRQ(ierr); 2398840c2264SHong Zhang ierr = VecAssemblyEnd(vo_nz);CHKERRQ(ierr); 2399840c2264SHong Zhang 2400840c2264SHong Zhang ierr = VecGetArray(vd_nz,&vdnz);CHKERRQ(ierr); 2401840c2264SHong Zhang ierr = VecGetArray(vo_nz,&vonz);CHKERRQ(ierr); 2402840c2264SHong Zhang for (j=0; j<localSize; j++) { 2403e102a522SHong Zhang dnnz[j] = (PetscInt)PetscRealPart(vdnz[j]); 2404e102a522SHong Zhang onnz[j] = (PetscInt)PetscRealPart(vonz[j]); 2405840c2264SHong Zhang } 2406840c2264SHong Zhang ierr = VecRestoreArray(vd_nz,&vdnz);CHKERRQ(ierr); 2407840c2264SHong Zhang ierr = VecRestoreArray(vo_nz,&vonz);CHKERRQ(ierr); 2408840c2264SHong Zhang ierr = VecDestroy(&vd_nz);CHKERRQ(ierr); 2409840c2264SHong Zhang ierr = VecDestroy(&vo_nz);CHKERRQ(ierr); 2410840c2264SHong Zhang 24115cf7da58SHong Zhang ierr = MatSeqAIJSetPreallocation(*J,0,dnnz);CHKERRQ(ierr); 24125cf7da58SHong Zhang ierr = MatMPIAIJSetPreallocation(*J,0,dnnz,0,onnz);CHKERRQ(ierr); 24135cf7da58SHong Zhang ierr = MatSetOption(*J,MAT_NEW_NONZERO_ALLOCATION_ERR,PETSC_TRUE);CHKERRQ(ierr); 24145cf7da58SHong Zhang 24155cf7da58SHong Zhang ierr = PetscFree2(dnnz,onnz);CHKERRQ(ierr); 24165cf7da58SHong Zhang 241789898e50SHong Zhang /* (2) Set matrix entries for edges */ 241889898e50SHong Zhang /*----------------------------------*/ 24191ad426b7SHong Zhang for (e=eStart; e<eEnd; e++) { 2420bfbc38dcSHong Zhang /* Get row indices */ 24212bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,e,ALL_COMPONENTS,&rstart);CHKERRQ(ierr); 24222bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,e,&nrows);CHKERRQ(ierr); 24234b976069SHong Zhang if (nrows) { 242417df6e9eSHong Zhang for (j=0; j<nrows; j++) rows[j] = j + rstart; 24251ad426b7SHong Zhang 2426a5b23f4aSJose E. Roman /* Set matrix entries for connected vertices */ 2427d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,e,&cone);CHKERRQ(ierr); 2428bfbc38dcSHong Zhang for (v=0; v<2; v++) { 24292bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,cone[v],ALL_COMPONENTS,&cstart);CHKERRQ(ierr); 24302bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,cone[v],&ncols);CHKERRQ(ierr); 24313e97b6e8SHong Zhang 24328675203cSHong Zhang if (network->Je) { 2433a4e85ca8SHong Zhang Juser = network->Je[3*e+1+v]; /* Jacobian(e,v) */ 24348675203cSHong Zhang } else Juser = NULL; 2435a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,ncols,cstart,J);CHKERRQ(ierr); 2436bfbc38dcSHong Zhang } 243717df6e9eSHong Zhang 2438bfbc38dcSHong Zhang /* Set matrix entries for edge self */ 24393e97b6e8SHong Zhang cstart = rstart; 24408675203cSHong Zhang if (network->Je) { 2441a4e85ca8SHong Zhang Juser = network->Je[3*e]; /* Jacobian(e,e) */ 24428675203cSHong Zhang } else Juser = NULL; 2443a4e85ca8SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows,nrows,cstart,J);CHKERRQ(ierr); 24441ad426b7SHong Zhang } 24454b976069SHong Zhang } 24461ad426b7SHong Zhang 2447bfbc38dcSHong Zhang /* Set matrix entries for vertices */ 244883b2e829SHong Zhang /*---------------------------------*/ 24491ad426b7SHong Zhang for (v=vStart; v<vEnd; v++) { 2450bfbc38dcSHong Zhang /* Get row indices */ 24512bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&rstart);CHKERRQ(ierr); 24522bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,v,&nrows);CHKERRQ(ierr); 24534b976069SHong Zhang if (!nrows) continue; 2454596e729fSHong Zhang 2455bdcb62a2SHong Zhang ierr = DMNetworkIsGhostVertex(dm,v,&ghost);CHKERRQ(ierr); 2456bdcb62a2SHong Zhang if (ghost) { 2457bdcb62a2SHong Zhang ierr = PetscMalloc1(nrows,&rows_v);CHKERRQ(ierr); 2458bdcb62a2SHong Zhang } else { 2459bdcb62a2SHong Zhang rows_v = rows; 2460bdcb62a2SHong Zhang } 2461bdcb62a2SHong Zhang for (j=0; j<nrows; j++) rows_v[j] = j + rstart; 2462596e729fSHong Zhang 2463bfbc38dcSHong Zhang /* Get supporting edges and connected vertices */ 2464596e729fSHong Zhang ierr = DMNetworkGetSupportingEdges(dm,v,&nedges,&edges);CHKERRQ(ierr); 2465596e729fSHong Zhang 2466596e729fSHong Zhang for (e=0; e<nedges; e++) { 2467bfbc38dcSHong Zhang /* Supporting edges */ 24682bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,edges[e],ALL_COMPONENTS,&cstart);CHKERRQ(ierr); 24692bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,edges[e],&ncols);CHKERRQ(ierr); 2470596e729fSHong Zhang 24718675203cSHong Zhang if (network->Jv) { 2472a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+1]; /* Jacobian(v,e) */ 24738675203cSHong Zhang } else Juser = NULL; 2474bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 2475596e729fSHong Zhang 2476bfbc38dcSHong Zhang /* Connected vertices */ 2477d842c372SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,edges[e],&cone);CHKERRQ(ierr); 24782a945128SHong Zhang vc = (v == cone[0]) ? cone[1]:cone[0]; 24792a945128SHong Zhang 24802bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,vc,ALL_COMPONENTS,&cstart);CHKERRQ(ierr); 24812bf73ac6SHong Zhang ierr = PetscSectionGetDof(network->DofSection,vc,&ncols);CHKERRQ(ierr); 2482a4e85ca8SHong Zhang 24838675203cSHong Zhang if (network->Jv) { 2484a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]+2*e+2]; /* Jacobian(v,vc) */ 24858675203cSHong Zhang } else Juser = NULL; 2486bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,ncols,cstart,J);CHKERRQ(ierr); 2487596e729fSHong Zhang } 2488596e729fSHong Zhang 2489bfbc38dcSHong Zhang /* Set matrix entries for vertex self */ 24901ad426b7SHong Zhang if (!ghost) { 24912bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,v,ALL_COMPONENTS,&cstart);CHKERRQ(ierr); 24928675203cSHong Zhang if (network->Jv) { 2493a4e85ca8SHong Zhang Juser = network->Jv[vptr[v-vStart]]; /* Jacobian(v,v) */ 24948675203cSHong Zhang } else Juser = NULL; 2495bdcb62a2SHong Zhang ierr = MatSetblock_private(Juser,nrows,rows_v,nrows,cstart,J);CHKERRQ(ierr); 2496bdcb62a2SHong Zhang } 2497bdcb62a2SHong Zhang if (ghost) { 2498bdcb62a2SHong Zhang ierr = PetscFree(rows_v);CHKERRQ(ierr); 2499bdcb62a2SHong Zhang } 25001ad426b7SHong Zhang } 2501a4e85ca8SHong Zhang ierr = PetscFree(rows);CHKERRQ(ierr); 2502bdcb62a2SHong Zhang 25031ad426b7SHong Zhang ierr = MatAssemblyBegin(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 25041ad426b7SHong Zhang ierr = MatAssemblyEnd(*J,MAT_FINAL_ASSEMBLY);CHKERRQ(ierr); 2505dd6f46cdSHong Zhang 25065f2c45f1SShri Abhyankar ierr = MatSetDM(*J,dm);CHKERRQ(ierr); 25075f2c45f1SShri Abhyankar PetscFunctionReturn(0); 25085f2c45f1SShri Abhyankar } 25095f2c45f1SShri Abhyankar 25105f2c45f1SShri Abhyankar PetscErrorCode DMDestroy_Network(DM dm) 25115f2c45f1SShri Abhyankar { 25125f2c45f1SShri Abhyankar PetscErrorCode ierr; 25135f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 251454dfd506SHong Zhang PetscInt j,np; 25155f2c45f1SShri Abhyankar 25165f2c45f1SShri Abhyankar PetscFunctionBegin; 25178415c774SShri Abhyankar if (--network->refct > 0) PetscFunctionReturn(0); 251883b2e829SHong Zhang if (network->Je) { 251983b2e829SHong Zhang ierr = PetscFree(network->Je);CHKERRQ(ierr); 252083b2e829SHong Zhang } 252183b2e829SHong Zhang if (network->Jv) { 2522883e35e8SHong Zhang ierr = PetscFree(network->Jvptr);CHKERRQ(ierr); 252383b2e829SHong Zhang ierr = PetscFree(network->Jv);CHKERRQ(ierr); 25241ad426b7SHong Zhang } 252513c2a604SAdrian Maldonado 252613c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->vertex.mapping);CHKERRQ(ierr); 252713c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.DofSection);CHKERRQ(ierr); 252813c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->vertex.GlobalDofSection);CHKERRQ(ierr); 2529f5427c60SHong Zhang if (network->vltog) { 2530f5427c60SHong Zhang ierr = PetscFree(network->vltog);CHKERRQ(ierr); 2531f5427c60SHong Zhang } 253213c2a604SAdrian Maldonado if (network->vertex.sf) { 253313c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->vertex.sf);CHKERRQ(ierr); 253413c2a604SAdrian Maldonado } 253513c2a604SAdrian Maldonado /* edge */ 253613c2a604SAdrian Maldonado ierr = ISLocalToGlobalMappingDestroy(&network->edge.mapping);CHKERRQ(ierr); 253713c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.DofSection);CHKERRQ(ierr); 253813c2a604SAdrian Maldonado ierr = PetscSectionDestroy(&network->edge.GlobalDofSection);CHKERRQ(ierr); 253913c2a604SAdrian Maldonado if (network->edge.sf) { 254013c2a604SAdrian Maldonado ierr = PetscSFDestroy(&network->edge.sf);CHKERRQ(ierr); 254113c2a604SAdrian Maldonado } 25425f2c45f1SShri Abhyankar ierr = DMDestroy(&network->plex);CHKERRQ(ierr); 25435f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DataSection);CHKERRQ(ierr); 25445f2c45f1SShri Abhyankar ierr = PetscSectionDestroy(&network->DofSection);CHKERRQ(ierr); 254583b2e829SHong Zhang 25462bf73ac6SHong Zhang for (j=0; j<network->Nsvtx; j++) { 25472bf73ac6SHong Zhang ierr = PetscFree(network->svtx[j].sv);CHKERRQ(ierr); 25482bf73ac6SHong Zhang } 25492bf73ac6SHong Zhang if (network->svtx) {ierr = PetscFree(network->svtx);CHKERRQ(ierr);} 2550f11a936eSBarry Smith ierr = PetscFree2(network->subnetedge,network->subnetvtx);CHKERRQ(ierr); 2551caf410d2SHong Zhang 25522bf73ac6SHong Zhang ierr = PetscTableDestroy(&network->svtable);CHKERRQ(ierr); 2553e2aaf10cSShri Abhyankar ierr = PetscFree(network->subnet);CHKERRQ(ierr); 255454dfd506SHong Zhang ierr = PetscFree(network->component);CHKERRQ(ierr); 25555f2c45f1SShri Abhyankar ierr = PetscFree(network->componentdataarray);CHKERRQ(ierr); 255654dfd506SHong Zhang 255754dfd506SHong Zhang if (network->header) { 255854dfd506SHong Zhang np = network->pEnd - network->pStart; 255954dfd506SHong Zhang for (j=0; j < np; j++) { 256054dfd506SHong Zhang ierr = PetscFree5(network->header[j].size,network->header[j].key,network->header[j].offset,network->header[j].nvar,network->header[j].offsetvarrel);CHKERRQ(ierr); 256154dfd506SHong Zhang ierr = PetscFree(network->cvalue[j].data);CHKERRQ(ierr); 256254dfd506SHong Zhang } 2563caf410d2SHong Zhang ierr = PetscFree2(network->header,network->cvalue);CHKERRQ(ierr); 256454dfd506SHong Zhang } 25655f2c45f1SShri Abhyankar ierr = PetscFree(network);CHKERRQ(ierr); 25665f2c45f1SShri Abhyankar PetscFunctionReturn(0); 25675f2c45f1SShri Abhyankar } 25685f2c45f1SShri Abhyankar 25695f2c45f1SShri Abhyankar PetscErrorCode DMView_Network(DM dm,PetscViewer viewer) 25705f2c45f1SShri Abhyankar { 25715f2c45f1SShri Abhyankar PetscErrorCode ierr; 2572caf410d2SHong Zhang PetscBool iascii; 2573caf410d2SHong Zhang PetscMPIInt rank; 25745f2c45f1SShri Abhyankar 25755f2c45f1SShri Abhyankar PetscFunctionBegin; 25769ace16cdSJacob Faibussowitsch PetscAssertFalse(!dm->setupcalled,PetscObjectComm((PetscObject)dm), PETSC_ERR_ARG_WRONGSTATE,"Must call DMSetUp() first"); 2577ffc4695bSBarry Smith ierr = MPI_Comm_rank(PetscObjectComm((PetscObject)dm),&rank);CHKERRMPI(ierr); 2578caf410d2SHong Zhang PetscValidHeaderSpecific(dm,DM_CLASSID, 1); 2579caf410d2SHong Zhang PetscValidHeaderSpecific(viewer, PETSC_VIEWER_CLASSID, 2); 2580caf410d2SHong Zhang ierr = PetscObjectTypeCompare((PetscObject)viewer,PETSCVIEWERASCII,&iascii);CHKERRQ(ierr); 2581caf410d2SHong Zhang if (iascii) { 2582caf410d2SHong Zhang const PetscInt *cone,*vtx,*edges; 25832bf73ac6SHong Zhang PetscInt vfrom,vto,i,j,nv,ne,ncv,p,nsubnet; 25842bf73ac6SHong Zhang DM_Network *network = (DM_Network*)dm->data; 2585caf410d2SHong Zhang 25862bf73ac6SHong Zhang nsubnet = network->Nsubnet; /* num of subnetworks */ 2587dd400576SPatrick Sanan if (rank == 0) { 25882bf73ac6SHong Zhang ierr = PetscPrintf(PETSC_COMM_SELF," NSubnets: %D; NEdges: %D; NVertices: %D; NSharedVertices: %D.\n",nsubnet,network->NEdges,network->NVertices,network->Nsvtx);CHKERRQ(ierr); 25892bf73ac6SHong Zhang } 25902bf73ac6SHong Zhang 25912bf73ac6SHong Zhang ierr = DMNetworkGetSharedVertices(dm,&ncv,&vtx);CHKERRQ(ierr); 2592caf410d2SHong Zhang ierr = PetscViewerASCIIPushSynchronized(viewer);CHKERRQ(ierr); 25932bf73ac6SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " [%d] nEdges: %D; nVertices: %D; nSharedVertices: %D\n",rank,network->nEdges,network->nVertices,ncv);CHKERRQ(ierr); 2594caf410d2SHong Zhang 2595caf410d2SHong Zhang for (i=0; i<nsubnet; i++) { 25962bf73ac6SHong Zhang ierr = DMNetworkGetSubnetwork(dm,i,&nv,&ne,&vtx,&edges);CHKERRQ(ierr); 2597caf410d2SHong Zhang if (ne) { 25982bf73ac6SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " Subnet %D: nEdges %D, nVertices(include shared vertices) %D\n",i,ne,nv);CHKERRQ(ierr); 2599caf410d2SHong Zhang for (j=0; j<ne; j++) { 2600caf410d2SHong Zhang p = edges[j]; 2601caf410d2SHong Zhang ierr = DMNetworkGetConnectedVertices(dm,p,&cone);CHKERRQ(ierr); 2602caf410d2SHong Zhang ierr = DMNetworkGetGlobalVertexIndex(dm,cone[0],&vfrom);CHKERRQ(ierr); 2603caf410d2SHong Zhang ierr = DMNetworkGetGlobalVertexIndex(dm,cone[1],&vto);CHKERRQ(ierr); 2604caf410d2SHong Zhang ierr = DMNetworkGetGlobalEdgeIndex(dm,edges[j],&p);CHKERRQ(ierr); 2605caf410d2SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " edge %D: %D ----> %D\n",p,vfrom,vto);CHKERRQ(ierr); 2606caf410d2SHong Zhang } 2607caf410d2SHong Zhang } 2608caf410d2SHong Zhang } 26092bf73ac6SHong Zhang 26102bf73ac6SHong Zhang /* Shared vertices */ 26112bf73ac6SHong Zhang ierr = DMNetworkGetSharedVertices(dm,&ncv,&vtx);CHKERRQ(ierr); 26122bf73ac6SHong Zhang if (ncv) { 26132bf73ac6SHong Zhang SVtx *svtx = network->svtx; 26142bf73ac6SHong Zhang PetscInt gidx,svtx_idx,nvto,vfrom_net,vfrom_idx,*svto; 26152bf73ac6SHong Zhang PetscBool ghost; 26162bf73ac6SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " SharedVertices:\n");CHKERRQ(ierr); 26172bf73ac6SHong Zhang for (i=0; i<ncv; i++) { 26182bf73ac6SHong Zhang ierr = DMNetworkIsGhostVertex(dm,vtx[i],&ghost);CHKERRQ(ierr); 26192bf73ac6SHong Zhang if (ghost) continue; 26202bf73ac6SHong Zhang 26212bf73ac6SHong Zhang ierr = DMNetworkGetGlobalVertexIndex(dm,vtx[i],&gidx);CHKERRQ(ierr); 26222bf73ac6SHong Zhang ierr = PetscTableFind(network->svtable,gidx+1,&svtx_idx);CHKERRQ(ierr); 26232bf73ac6SHong Zhang svtx_idx--; 26242bf73ac6SHong Zhang nvto = svtx[svtx_idx].n; 26252bf73ac6SHong Zhang 26262bf73ac6SHong Zhang vfrom_net = svtx[svtx_idx].sv[0]; 26272bf73ac6SHong Zhang vfrom_idx = svtx[svtx_idx].sv[1]; 26282bf73ac6SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " svtx %D: global index %D, subnet[%D].%D ---->\n",i,gidx,vfrom_net,vfrom_idx);CHKERRQ(ierr); 26292bf73ac6SHong Zhang for (j=1; j<nvto; j++) { 26302bf73ac6SHong Zhang svto = svtx[svtx_idx].sv + 2*j; 26312bf73ac6SHong Zhang ierr = PetscViewerASCIISynchronizedPrintf(viewer, " ----> subnet[%D].%D\n",svto[0],svto[1]);CHKERRQ(ierr); 2632caf410d2SHong Zhang } 2633caf410d2SHong Zhang } 2634caf410d2SHong Zhang } 2635caf410d2SHong Zhang ierr = PetscViewerFlush(viewer);CHKERRQ(ierr); 2636caf410d2SHong Zhang ierr = PetscViewerASCIIPopSynchronized(viewer);CHKERRQ(ierr); 263798921bdaSJacob Faibussowitsch } else SETERRQ(PetscObjectComm((PetscObject) dm), PETSC_ERR_SUP, "Viewer type %s not yet supported for DMNetwork writing", ((PetscObject)viewer)->type_name); 26385f2c45f1SShri Abhyankar PetscFunctionReturn(0); 26395f2c45f1SShri Abhyankar } 26405f2c45f1SShri Abhyankar 26415f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalBegin_Network(DM dm, Vec g, InsertMode mode, Vec l) 26425f2c45f1SShri Abhyankar { 26435f2c45f1SShri Abhyankar PetscErrorCode ierr; 26445f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 26455f2c45f1SShri Abhyankar 26465f2c45f1SShri Abhyankar PetscFunctionBegin; 26475f2c45f1SShri Abhyankar ierr = DMGlobalToLocalBegin(network->plex,g,mode,l);CHKERRQ(ierr); 26485f2c45f1SShri Abhyankar PetscFunctionReturn(0); 26495f2c45f1SShri Abhyankar } 26505f2c45f1SShri Abhyankar 26515f2c45f1SShri Abhyankar PetscErrorCode DMGlobalToLocalEnd_Network(DM dm, Vec g, InsertMode mode, Vec l) 26525f2c45f1SShri Abhyankar { 26535f2c45f1SShri Abhyankar PetscErrorCode ierr; 26545f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 26555f2c45f1SShri Abhyankar 26565f2c45f1SShri Abhyankar PetscFunctionBegin; 26575f2c45f1SShri Abhyankar ierr = DMGlobalToLocalEnd(network->plex,g,mode,l);CHKERRQ(ierr); 26585f2c45f1SShri Abhyankar PetscFunctionReturn(0); 26595f2c45f1SShri Abhyankar } 26605f2c45f1SShri Abhyankar 26615f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalBegin_Network(DM dm, Vec l, InsertMode mode, Vec g) 26625f2c45f1SShri Abhyankar { 26635f2c45f1SShri Abhyankar PetscErrorCode ierr; 26645f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 26655f2c45f1SShri Abhyankar 26665f2c45f1SShri Abhyankar PetscFunctionBegin; 26675f2c45f1SShri Abhyankar ierr = DMLocalToGlobalBegin(network->plex,l,mode,g);CHKERRQ(ierr); 26685f2c45f1SShri Abhyankar PetscFunctionReturn(0); 26695f2c45f1SShri Abhyankar } 26705f2c45f1SShri Abhyankar 26715f2c45f1SShri Abhyankar PetscErrorCode DMLocalToGlobalEnd_Network(DM dm, Vec l, InsertMode mode, Vec g) 26725f2c45f1SShri Abhyankar { 26735f2c45f1SShri Abhyankar PetscErrorCode ierr; 26745f2c45f1SShri Abhyankar DM_Network *network = (DM_Network*)dm->data; 26755f2c45f1SShri Abhyankar 26765f2c45f1SShri Abhyankar PetscFunctionBegin; 26775f2c45f1SShri Abhyankar ierr = DMLocalToGlobalEnd(network->plex,l,mode,g);CHKERRQ(ierr); 26785f2c45f1SShri Abhyankar PetscFunctionReturn(0); 26795f2c45f1SShri Abhyankar } 268022bbedd7SHong Zhang 268122bbedd7SHong Zhang /*@ 268264238763SRodolfo Oliveira DMNetworkGetVertexLocalToGlobalOrdering - Get vertex global index 268322bbedd7SHong Zhang 268422bbedd7SHong Zhang Not collective 268522bbedd7SHong Zhang 26867a7aea1fSJed Brown Input Parameters: 268722bbedd7SHong Zhang + dm - the dm object 268822bbedd7SHong Zhang - vloc - local vertex ordering, start from 0 268922bbedd7SHong Zhang 26907a7aea1fSJed Brown Output Parameters: 2691f0fc11ceSJed Brown . vg - global vertex ordering, start from 0 269222bbedd7SHong Zhang 269397bb938eSShri Abhyankar Level: advanced 269422bbedd7SHong Zhang 269522bbedd7SHong Zhang .seealso: DMNetworkSetVertexLocalToGlobalOrdering() 269622bbedd7SHong Zhang @*/ 269722bbedd7SHong Zhang PetscErrorCode DMNetworkGetVertexLocalToGlobalOrdering(DM dm,PetscInt vloc,PetscInt *vg) 269822bbedd7SHong Zhang { 269922bbedd7SHong Zhang DM_Network *network = (DM_Network*)dm->data; 270022bbedd7SHong Zhang PetscInt *vltog = network->vltog; 270122bbedd7SHong Zhang 270222bbedd7SHong Zhang PetscFunctionBegin; 27039ace16cdSJacob Faibussowitsch PetscAssertFalse(!vltog,PetscObjectComm((PetscObject)dm),PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkSetVertexLocalToGlobalOrdering() first"); 270422bbedd7SHong Zhang *vg = vltog[vloc]; 270522bbedd7SHong Zhang PetscFunctionReturn(0); 270622bbedd7SHong Zhang } 270722bbedd7SHong Zhang 270822bbedd7SHong Zhang /*@ 270964238763SRodolfo Oliveira DMNetworkSetVertexLocalToGlobalOrdering - Create and setup vertex local to global map 271022bbedd7SHong Zhang 271122bbedd7SHong Zhang Collective 271222bbedd7SHong Zhang 271322bbedd7SHong Zhang Input Parameters: 2714f0fc11ceSJed Brown . dm - the dm object 271522bbedd7SHong Zhang 271697bb938eSShri Abhyankar Level: advanced 271722bbedd7SHong Zhang 271863029d29SHong Zhang .seealso: DMNetworkGetGlobalVertexIndex() 271922bbedd7SHong Zhang @*/ 272022bbedd7SHong Zhang PetscErrorCode DMNetworkSetVertexLocalToGlobalOrdering(DM dm) 272122bbedd7SHong Zhang { 272222bbedd7SHong Zhang PetscErrorCode ierr; 272322bbedd7SHong Zhang DM_Network *network=(DM_Network*)dm->data; 272422bbedd7SHong Zhang MPI_Comm comm; 27252bf73ac6SHong Zhang PetscMPIInt rank,size,*displs=NULL,*recvcounts=NULL,remoterank; 272622bbedd7SHong Zhang PetscBool ghost; 272763029d29SHong Zhang PetscInt *vltog,nroots,nleaves,i,*vrange,k,N,lidx; 272822bbedd7SHong Zhang const PetscSFNode *iremote; 272922bbedd7SHong Zhang PetscSF vsf; 273063029d29SHong Zhang Vec Vleaves,Vleaves_seq; 273163029d29SHong Zhang VecScatter ctx; 273263029d29SHong Zhang PetscScalar *varr,val; 273363029d29SHong Zhang const PetscScalar *varr_read; 273422bbedd7SHong Zhang 273522bbedd7SHong Zhang PetscFunctionBegin; 273622bbedd7SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 2737ffc4695bSBarry Smith ierr = MPI_Comm_size(comm,&size);CHKERRMPI(ierr); 2738ffc4695bSBarry Smith ierr = MPI_Comm_rank(comm,&rank);CHKERRMPI(ierr); 273922bbedd7SHong Zhang 274022bbedd7SHong Zhang if (size == 1) { 274122bbedd7SHong Zhang nroots = network->vEnd - network->vStart; 274222bbedd7SHong Zhang ierr = PetscMalloc1(nroots, &vltog);CHKERRQ(ierr); 274322bbedd7SHong Zhang for (i=0; i<nroots; i++) vltog[i] = i; 274422bbedd7SHong Zhang network->vltog = vltog; 274522bbedd7SHong Zhang PetscFunctionReturn(0); 274622bbedd7SHong Zhang } 274722bbedd7SHong Zhang 27489ace16cdSJacob Faibussowitsch PetscAssertFalse(!network->distributecalled,comm, PETSC_ERR_ARG_WRONGSTATE,"Must call DMNetworkDistribute() first"); 274922bbedd7SHong Zhang if (network->vltog) { 275022bbedd7SHong Zhang ierr = PetscFree(network->vltog);CHKERRQ(ierr); 275122bbedd7SHong Zhang } 275222bbedd7SHong Zhang 275322bbedd7SHong Zhang ierr = DMNetworkSetSubMap_private(network->vStart,network->vEnd,&network->vertex.mapping);CHKERRQ(ierr); 275422bbedd7SHong Zhang ierr = PetscSFGetSubSF(network->plex->sf, network->vertex.mapping, &network->vertex.sf);CHKERRQ(ierr); 275522bbedd7SHong Zhang vsf = network->vertex.sf; 275622bbedd7SHong Zhang 27572bf73ac6SHong Zhang ierr = PetscMalloc3(size+1,&vrange,size,&displs,size,&recvcounts);CHKERRQ(ierr); 2758f5427c60SHong Zhang ierr = PetscSFGetGraph(vsf,&nroots,&nleaves,NULL,&iremote);CHKERRQ(ierr); 275922bbedd7SHong Zhang 276022bbedd7SHong Zhang for (i=0; i<size; i++) { displs[i] = i; recvcounts[i] = 1;} 276122bbedd7SHong Zhang 276222bbedd7SHong Zhang i = nroots - nleaves; /* local number of vertices, excluding ghosts */ 276322bbedd7SHong Zhang vrange[0] = 0; 2764ffc4695bSBarry Smith ierr = MPI_Allgatherv(&i,1,MPIU_INT,vrange+1,recvcounts,displs,MPIU_INT,comm);CHKERRMPI(ierr); 276567dd800bSHong Zhang for (i=2; i<size+1; i++) {vrange[i] += vrange[i-1];} 276622bbedd7SHong Zhang 276722bbedd7SHong Zhang ierr = PetscMalloc1(nroots, &vltog);CHKERRQ(ierr); 276822bbedd7SHong Zhang network->vltog = vltog; 276922bbedd7SHong Zhang 277063029d29SHong Zhang /* Set vltog for non-ghost vertices */ 277163029d29SHong Zhang k = 0; 277222bbedd7SHong Zhang for (i=0; i<nroots; i++) { 277322bbedd7SHong Zhang ierr = DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost);CHKERRQ(ierr); 277463029d29SHong Zhang if (ghost) continue; 277563029d29SHong Zhang vltog[i] = vrange[rank] + k++; 277622bbedd7SHong Zhang } 2777f5427c60SHong Zhang ierr = PetscFree3(vrange,displs,recvcounts);CHKERRQ(ierr); 277863029d29SHong Zhang 277963029d29SHong Zhang /* Set vltog for ghost vertices */ 278063029d29SHong Zhang /* (a) create parallel Vleaves and sequential Vleaves_seq to convert local iremote[*].index to global index */ 278163029d29SHong Zhang ierr = VecCreate(comm,&Vleaves);CHKERRQ(ierr); 278263029d29SHong Zhang ierr = VecSetSizes(Vleaves,2*nleaves,PETSC_DETERMINE);CHKERRQ(ierr); 278363029d29SHong Zhang ierr = VecSetFromOptions(Vleaves);CHKERRQ(ierr); 278463029d29SHong Zhang ierr = VecGetArray(Vleaves,&varr);CHKERRQ(ierr); 278563029d29SHong Zhang for (i=0; i<nleaves; i++) { 278663029d29SHong Zhang varr[2*i] = (PetscScalar)(iremote[i].rank); /* rank of remote process */ 278763029d29SHong Zhang varr[2*i+1] = (PetscScalar)(iremote[i].index); /* local index in remote process */ 278863029d29SHong Zhang } 278963029d29SHong Zhang ierr = VecRestoreArray(Vleaves,&varr);CHKERRQ(ierr); 279063029d29SHong Zhang 279163029d29SHong Zhang /* (b) scatter local info to remote processes via VecScatter() */ 279263029d29SHong Zhang ierr = VecScatterCreateToAll(Vleaves,&ctx,&Vleaves_seq);CHKERRQ(ierr); 279363029d29SHong Zhang ierr = VecScatterBegin(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 279463029d29SHong Zhang ierr = VecScatterEnd(ctx,Vleaves,Vleaves_seq,INSERT_VALUES,SCATTER_FORWARD);CHKERRQ(ierr); 279563029d29SHong Zhang 279663029d29SHong Zhang /* (c) convert local indices to global indices in parallel vector Vleaves */ 279763029d29SHong Zhang ierr = VecGetSize(Vleaves_seq,&N);CHKERRQ(ierr); 279863029d29SHong Zhang ierr = VecGetArrayRead(Vleaves_seq,&varr_read);CHKERRQ(ierr); 279963029d29SHong Zhang for (i=0; i<N; i+=2) { 28002e4cff2eSHong Zhang remoterank = (PetscMPIInt)PetscRealPart(varr_read[i]); 280163029d29SHong Zhang if (remoterank == rank) { 280263029d29SHong Zhang k = i+1; /* row number */ 28032e4cff2eSHong Zhang lidx = (PetscInt)PetscRealPart(varr_read[i+1]); 280463029d29SHong Zhang val = (PetscScalar)vltog[lidx]; /* global index for non-ghost vertex computed above */ 280563029d29SHong Zhang ierr = VecSetValues(Vleaves,1,&k,&val,INSERT_VALUES);CHKERRQ(ierr); 280663029d29SHong Zhang } 280763029d29SHong Zhang } 280863029d29SHong Zhang ierr = VecRestoreArrayRead(Vleaves_seq,&varr_read);CHKERRQ(ierr); 280963029d29SHong Zhang ierr = VecAssemblyBegin(Vleaves);CHKERRQ(ierr); 281063029d29SHong Zhang ierr = VecAssemblyEnd(Vleaves);CHKERRQ(ierr); 281163029d29SHong Zhang 281263029d29SHong Zhang /* (d) Set vltog for ghost vertices by copying local values of Vleaves */ 281363029d29SHong Zhang ierr = VecGetArrayRead(Vleaves,&varr_read);CHKERRQ(ierr); 281463029d29SHong Zhang k = 0; 281563029d29SHong Zhang for (i=0; i<nroots; i++) { 281663029d29SHong Zhang ierr = DMNetworkIsGhostVertex(dm,i+network->vStart,&ghost);CHKERRQ(ierr); 281763029d29SHong Zhang if (!ghost) continue; 28182e4cff2eSHong Zhang vltog[i] = (PetscInt)PetscRealPart(varr_read[2*k+1]); k++; 281963029d29SHong Zhang } 282063029d29SHong Zhang ierr = VecRestoreArrayRead(Vleaves,&varr_read);CHKERRQ(ierr); 282163029d29SHong Zhang 282263029d29SHong Zhang ierr = VecDestroy(&Vleaves);CHKERRQ(ierr); 282363029d29SHong Zhang ierr = VecDestroy(&Vleaves_seq);CHKERRQ(ierr); 282463029d29SHong Zhang ierr = VecScatterDestroy(&ctx);CHKERRQ(ierr); 282522bbedd7SHong Zhang PetscFunctionReturn(0); 282622bbedd7SHong Zhang } 282742dc13f1SHong Zhang 2828*9fbee547SJacob Faibussowitsch static inline PetscErrorCode DMISAddSize_private(DM_Network *network,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *nidx) 282942dc13f1SHong Zhang { 283042dc13f1SHong Zhang PetscErrorCode ierr; 283142dc13f1SHong Zhang PetscInt i,j,ncomps,nvar,key,offset=0; 283242dc13f1SHong Zhang DMNetworkComponentHeader header; 283342dc13f1SHong Zhang 283442dc13f1SHong Zhang PetscFunctionBegin; 283542dc13f1SHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 283642dc13f1SHong Zhang ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 283742dc13f1SHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offset); 283842dc13f1SHong Zhang 283942dc13f1SHong Zhang for (i=0; i<ncomps; i++) { 284042dc13f1SHong Zhang key = header->key[i]; 284142dc13f1SHong Zhang nvar = header->nvar[i]; 284242dc13f1SHong Zhang for (j=0; j<numkeys; j++) { 284342dc13f1SHong Zhang if (key == keys[j]) { 284442dc13f1SHong Zhang if (!blocksize || blocksize[j] == -1) { 284542dc13f1SHong Zhang *nidx += nvar; 284642dc13f1SHong Zhang } else { 284742dc13f1SHong Zhang *nidx += nselectedvar[j]*nvar/blocksize[j]; 284842dc13f1SHong Zhang } 284942dc13f1SHong Zhang } 285042dc13f1SHong Zhang } 285142dc13f1SHong Zhang } 285242dc13f1SHong Zhang PetscFunctionReturn(0); 285342dc13f1SHong Zhang } 285442dc13f1SHong Zhang 2855*9fbee547SJacob Faibussowitsch static inline PetscErrorCode DMISComputeIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx) 285642dc13f1SHong Zhang { 285742dc13f1SHong Zhang PetscErrorCode ierr; 285842dc13f1SHong Zhang PetscInt i,j,ncomps,nvar,key,offsetg,k,k1,offset=0; 285942dc13f1SHong Zhang DM_Network *network = (DM_Network*)dm->data; 286042dc13f1SHong Zhang DMNetworkComponentHeader header; 286142dc13f1SHong Zhang 286242dc13f1SHong Zhang PetscFunctionBegin; 286342dc13f1SHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 286442dc13f1SHong Zhang ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 286542dc13f1SHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offset); 286642dc13f1SHong Zhang 286742dc13f1SHong Zhang for (i=0; i<ncomps; i++) { 286842dc13f1SHong Zhang key = header->key[i]; 286942dc13f1SHong Zhang nvar = header->nvar[i]; 287042dc13f1SHong Zhang for (j=0; j<numkeys; j++) { 287142dc13f1SHong Zhang if (key != keys[j]) continue; 287242dc13f1SHong Zhang 287342dc13f1SHong Zhang ierr = DMNetworkGetGlobalVecOffset(dm,p,i,&offsetg);CHKERRQ(ierr); 287442dc13f1SHong Zhang if (!blocksize || blocksize[j] == -1) { 287542dc13f1SHong Zhang for (k=0; k<nvar; k++) idx[(*ii)++] = offsetg + k; 287642dc13f1SHong Zhang } else { 287742dc13f1SHong Zhang for (k=0; k<nvar; k+=blocksize[j]) { 287842dc13f1SHong Zhang for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetg + k + selectedvar[j][k1]; 287942dc13f1SHong Zhang } 288042dc13f1SHong Zhang } 288142dc13f1SHong Zhang } 288242dc13f1SHong Zhang } 288342dc13f1SHong Zhang PetscFunctionReturn(0); 288442dc13f1SHong Zhang } 288542dc13f1SHong Zhang 288642dc13f1SHong Zhang /*@ 288742dc13f1SHong Zhang DMNetworkCreateIS - Create an index set object from the global vector of the network 288842dc13f1SHong Zhang 288942dc13f1SHong Zhang Collective 289042dc13f1SHong Zhang 289142dc13f1SHong Zhang Input Parameters: 289242dc13f1SHong Zhang + dm - DMNetwork object 289342dc13f1SHong Zhang . numkeys - number of keys 289442dc13f1SHong Zhang . keys - array of keys that define the components of the variables you wish to extract 289542dc13f1SHong Zhang . blocksize - block size of the variables associated to the component 289642dc13f1SHong Zhang . nselectedvar - number of variables in each block to select 289742dc13f1SHong Zhang - selectedvar - the offset into the block of each variable in each block to select 289842dc13f1SHong Zhang 289942dc13f1SHong Zhang Output Parameters: 290042dc13f1SHong Zhang . is - the index set 290142dc13f1SHong Zhang 290242dc13f1SHong Zhang Level: Advanced 290342dc13f1SHong Zhang 290442dc13f1SHong Zhang Notes: 290542dc13f1SHong Zhang Use blocksize[i] of -1 to indicate select all the variables of the i-th component, for which nselectvar[i] and selectedvar[i] are ignored. Use NULL, NULL, NULL to indicate for all selected components one wishes to obtain all the values of that component. For example, DMNetworkCreateIS(dm,1,&key,NULL,NULL,NULL,&is) will return an is that extracts all the variables for the 0-th component. 290642dc13f1SHong Zhang 290742dc13f1SHong Zhang .seealso: DMNetworkCreate(), ISCreateGeneral(), DMNetworkCreateLocalIS() 290842dc13f1SHong Zhang @*/ 290942dc13f1SHong Zhang PetscErrorCode DMNetworkCreateIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is) 291042dc13f1SHong Zhang { 291142dc13f1SHong Zhang PetscErrorCode ierr; 291242dc13f1SHong Zhang MPI_Comm comm; 291342dc13f1SHong Zhang DM_Network *network = (DM_Network*)dm->data; 291442dc13f1SHong Zhang PetscInt i,p,estart,eend,vstart,vend,nidx,*idx; 291542dc13f1SHong Zhang PetscBool ghost; 291642dc13f1SHong Zhang 291742dc13f1SHong Zhang PetscFunctionBegin; 291842dc13f1SHong Zhang ierr = PetscObjectGetComm((PetscObject)dm,&comm);CHKERRQ(ierr); 291942dc13f1SHong Zhang 292042dc13f1SHong Zhang /* Check input parameters */ 292142dc13f1SHong Zhang for (i=0; i<numkeys; i++) { 292242dc13f1SHong Zhang if (!blocksize || blocksize[i] == -1) continue; 29239ace16cdSJacob Faibussowitsch PetscAssertFalse(nselectedvar[i] > blocksize[i],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"number of selectedvariables %D cannot be larger than blocksize %D",nselectedvar[i],blocksize[i]); 292442dc13f1SHong Zhang } 292542dc13f1SHong Zhang 292642dc13f1SHong Zhang ierr = DMNetworkGetEdgeRange(dm,&estart,&eend);CHKERRQ(ierr); 292742dc13f1SHong Zhang ierr = DMNetworkGetVertexRange(dm,&vstart,&vend);CHKERRQ(ierr); 292842dc13f1SHong Zhang 292942dc13f1SHong Zhang /* Get local number of idx */ 293042dc13f1SHong Zhang nidx = 0; 293142dc13f1SHong Zhang for (p=estart; p<eend; p++) { 293242dc13f1SHong Zhang ierr = DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);CHKERRQ(ierr); 293342dc13f1SHong Zhang } 293442dc13f1SHong Zhang for (p=vstart; p<vend; p++) { 293542dc13f1SHong Zhang ierr = DMNetworkIsGhostVertex(dm,p,&ghost);CHKERRQ(ierr); 293642dc13f1SHong Zhang if (ghost) continue; 293742dc13f1SHong Zhang ierr = DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);CHKERRQ(ierr); 293842dc13f1SHong Zhang } 293942dc13f1SHong Zhang 294042dc13f1SHong Zhang /* Compute idx */ 294142dc13f1SHong Zhang ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 294242dc13f1SHong Zhang i = 0; 294342dc13f1SHong Zhang for (p=estart; p<eend; p++) { 294442dc13f1SHong Zhang ierr = DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);CHKERRQ(ierr); 294542dc13f1SHong Zhang } 294642dc13f1SHong Zhang for (p=vstart; p<vend; p++) { 294742dc13f1SHong Zhang ierr = DMNetworkIsGhostVertex(dm,p,&ghost);CHKERRQ(ierr); 294842dc13f1SHong Zhang if (ghost) continue; 294942dc13f1SHong Zhang ierr = DMISComputeIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);CHKERRQ(ierr); 295042dc13f1SHong Zhang } 295142dc13f1SHong Zhang 295242dc13f1SHong Zhang /* Create is */ 295342dc13f1SHong Zhang ierr = ISCreateGeneral(comm,nidx,idx,PETSC_COPY_VALUES,is);CHKERRQ(ierr); 295442dc13f1SHong Zhang ierr = PetscFree(idx);CHKERRQ(ierr); 295542dc13f1SHong Zhang PetscFunctionReturn(0); 295642dc13f1SHong Zhang } 295742dc13f1SHong Zhang 2958*9fbee547SJacob Faibussowitsch static inline PetscErrorCode DMISComputeLocalIdx_private(DM dm,PetscInt p,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],PetscInt *ii,PetscInt *idx) 295942dc13f1SHong Zhang { 296042dc13f1SHong Zhang PetscErrorCode ierr; 296142dc13f1SHong Zhang PetscInt i,j,ncomps,nvar,key,offsetl,k,k1,offset=0; 296242dc13f1SHong Zhang DM_Network *network = (DM_Network*)dm->data; 296342dc13f1SHong Zhang DMNetworkComponentHeader header; 296442dc13f1SHong Zhang 296542dc13f1SHong Zhang PetscFunctionBegin; 296642dc13f1SHong Zhang ierr = PetscSectionGetOffset(network->DataSection,p,&offset);CHKERRQ(ierr); 296742dc13f1SHong Zhang ncomps = ((DMNetworkComponentHeader)(network->componentdataarray+offset))->ndata; 296842dc13f1SHong Zhang header = (DMNetworkComponentHeader)(network->componentdataarray+offset); 296942dc13f1SHong Zhang 297042dc13f1SHong Zhang for (i=0; i<ncomps; i++) { 297142dc13f1SHong Zhang key = header->key[i]; 297242dc13f1SHong Zhang nvar = header->nvar[i]; 297342dc13f1SHong Zhang for (j=0; j<numkeys; j++) { 297442dc13f1SHong Zhang if (key != keys[j]) continue; 297542dc13f1SHong Zhang 297642dc13f1SHong Zhang ierr = DMNetworkGetLocalVecOffset(dm,p,i,&offsetl);CHKERRQ(ierr); 297742dc13f1SHong Zhang if (!blocksize || blocksize[j] == -1) { 297842dc13f1SHong Zhang for (k=0; k<nvar; k++) idx[(*ii)++] = offsetl + k; 297942dc13f1SHong Zhang } else { 298042dc13f1SHong Zhang for (k=0; k<nvar; k+=blocksize[j]) { 298142dc13f1SHong Zhang for (k1=0; k1<nselectedvar[j]; k1++) idx[(*ii)++] = offsetl + k + selectedvar[j][k1]; 298242dc13f1SHong Zhang } 298342dc13f1SHong Zhang } 298442dc13f1SHong Zhang } 298542dc13f1SHong Zhang } 298642dc13f1SHong Zhang PetscFunctionReturn(0); 298742dc13f1SHong Zhang } 298842dc13f1SHong Zhang 298942dc13f1SHong Zhang /*@ 299042dc13f1SHong Zhang DMNetworkCreateLocalIS - Create an index set object from the local vector of the network 299142dc13f1SHong Zhang 299242dc13f1SHong Zhang Not collective 299342dc13f1SHong Zhang 299442dc13f1SHong Zhang Input Parameters: 299542dc13f1SHong Zhang + dm - DMNetwork object 299642dc13f1SHong Zhang . numkeys - number of keys 299742dc13f1SHong Zhang . keys - array of keys that define the components of the variables you wish to extract 299842dc13f1SHong Zhang . blocksize - block size of the variables associated to the component 299942dc13f1SHong Zhang . nselectedvar - number of variables in each block to select 300042dc13f1SHong Zhang - selectedvar - the offset into the block of each variable in each block to select 300142dc13f1SHong Zhang 300242dc13f1SHong Zhang Output Parameters: 300342dc13f1SHong Zhang . is - the index set 300442dc13f1SHong Zhang 300542dc13f1SHong Zhang Level: Advanced 300642dc13f1SHong Zhang 300742dc13f1SHong Zhang Notes: 300842dc13f1SHong Zhang Use blocksize[i] of -1 to indicate select all the variables of the i-th component, for which nselectvar[i] and selectedvar[i] are ignored. Use NULL, NULL, NULL to indicate for all selected components one wishes to obtain all the values of that component. For example, DMNetworkCreateLocalIS(dm,1,&key,NULL,NULL,NULL,&is) will return an is that extracts all the variables for the 0-th component. 300942dc13f1SHong Zhang 301042dc13f1SHong Zhang .seealso: DMNetworkCreate(), DMNetworkCreateIS, ISCreateGeneral() 301142dc13f1SHong Zhang @*/ 301242dc13f1SHong Zhang PetscErrorCode DMNetworkCreateLocalIS(DM dm,PetscInt numkeys,PetscInt keys[],PetscInt blocksize[],PetscInt nselectedvar[],PetscInt *selectedvar[],IS *is) 301342dc13f1SHong Zhang { 301442dc13f1SHong Zhang PetscErrorCode ierr; 301542dc13f1SHong Zhang DM_Network *network = (DM_Network*)dm->data; 301642dc13f1SHong Zhang PetscInt i,p,pstart,pend,nidx,*idx; 301742dc13f1SHong Zhang 301842dc13f1SHong Zhang PetscFunctionBegin; 301942dc13f1SHong Zhang /* Check input parameters */ 302042dc13f1SHong Zhang for (i=0; i<numkeys; i++) { 302142dc13f1SHong Zhang if (!blocksize || blocksize[i] == -1) continue; 30229ace16cdSJacob Faibussowitsch PetscAssertFalse(nselectedvar[i] > blocksize[i],PETSC_COMM_SELF,PETSC_ERR_ARG_OUTOFRANGE,"number of selectedvariables %D cannot be larger than blocksize %D",nselectedvar[i],blocksize[i]); 302342dc13f1SHong Zhang } 302442dc13f1SHong Zhang 302542dc13f1SHong Zhang pstart = network->pStart; 302642dc13f1SHong Zhang pend = network->pEnd; 302742dc13f1SHong Zhang 302842dc13f1SHong Zhang /* Get local number of idx */ 302942dc13f1SHong Zhang nidx = 0; 303042dc13f1SHong Zhang for (p=pstart; p<pend; p++) { 303142dc13f1SHong Zhang ierr = DMISAddSize_private(network,p,numkeys,keys,blocksize,nselectedvar,&nidx);CHKERRQ(ierr); 303242dc13f1SHong Zhang } 303342dc13f1SHong Zhang 303442dc13f1SHong Zhang /* Compute local idx */ 303542dc13f1SHong Zhang ierr = PetscMalloc1(nidx,&idx);CHKERRQ(ierr); 303642dc13f1SHong Zhang i = 0; 303742dc13f1SHong Zhang for (p=pstart; p<pend; p++) { 303842dc13f1SHong Zhang ierr = DMISComputeLocalIdx_private(dm,p,numkeys,keys,blocksize,nselectedvar,selectedvar,&i,idx);CHKERRQ(ierr); 303942dc13f1SHong Zhang } 304042dc13f1SHong Zhang 304142dc13f1SHong Zhang /* Create is */ 304242dc13f1SHong Zhang ierr = ISCreateGeneral(PETSC_COMM_SELF,nidx,idx,PETSC_COPY_VALUES,is);CHKERRQ(ierr); 304342dc13f1SHong Zhang ierr = PetscFree(idx);CHKERRQ(ierr); 304442dc13f1SHong Zhang PetscFunctionReturn(0); 304542dc13f1SHong Zhang } 3046