#ifdef PARALLEL #include #include "paralleltop.hpp" namespace ngcore { template auto Max (FlatArray array) -> T { T max = std::numeric_limits::min(); for (auto & v : array) if (v > max) max = v; return max; } /* template auto Max (FlatArray array, TB initial) -> T { T max = initial; for (auto & v : array) if (v > max) max = v; return max; } */ } namespace netgen { ParallelMeshTopology :: ParallelMeshTopology (const Mesh & amesh) : mesh(amesh) { is_updated = false; } ParallelMeshTopology :: ~ParallelMeshTopology () { ; } void ParallelMeshTopology :: Reset () { *testout << "ParallelMeshTopology::Reset" << endl; if ( mesh.GetCommunicator().Size() == 1 ) return; int ned = mesh.GetTopology().GetNEdges(); int nfa = mesh.GetTopology().GetNFaces(); if (glob_edge.Size() != ned) { glob_edge.SetSize(ned); glob_face.SetSize(nfa); glob_edge = -1; glob_face = -1; loc2distedge.ChangeSize (ned); loc2distface.ChangeSize (nfa); } if (glob_vert.Size() != mesh.GetNV()) { SetNV(mesh.GetNV()); SetNE(mesh.GetNE()); } } void ParallelMeshTopology :: Print() const { ; } void ParallelMeshTopology :: EnumeratePointsGlobally () { auto comm = mesh.GetCommunicator(); auto rank = comm.Rank(); size_t oldnv = glob_vert.Size(); size_t nv = loc2distvert.Size(); *testout << "enumerate globally, loc2distvert.size = " << loc2distvert.Size() << ", glob_vert.size = " << glob_vert.Size() << endl; // *testout << "old glob_vert = " << endl << glob_vert << endl; if (rank == 0) nv = 0; IntRange newvr(oldnv, nv); // new vertex range glob_vert.SetSize (nv); glob_vert.Range(newvr) = -1; int num_master_points = 0; for (auto i : newvr) { auto dps = GetDistantPNums(i); // check sorted: for (int j = 0; j+1 < dps.Size(); j++) if (dps[j+1] < dps[j]) cout << "wrong sort" << endl; if (dps.Size() == 0 || dps[0] > comm.Rank()) glob_vert[i] = num_master_points++; } *testout << "nummaster = " << num_master_points << endl; Array first_master_point(comm.Size()); comm.AllGather (num_master_points, first_master_point); auto max_oldv = comm.AllReduce (Max (glob_vert.Range(0, oldnv)), MPI_MAX); if (comm.AllReduce (oldnv, MPI_SUM) == 0) max_oldv = PointIndex::BASE-1; size_t num_glob_points = max_oldv+1; // PointIndex::BASE; for (int i = 0; i < comm.Size(); i++) { int cur = first_master_point[i]; first_master_point[i] = num_glob_points; num_glob_points += cur; } for (auto i : newvr) if (glob_vert[i] != -1) glob_vert[i] += first_master_point[comm.Rank()]; // ScatterDofData (global_nums); Array nsend(comm.Size()), nrecv(comm.Size()); nsend = 0; nrecv = 0; /** Count send/recv size **/ for (auto i : newvr) { auto dps = GetDistantPNums(i); if (!dps.Size()) continue; if (rank < dps[0]) for(auto p:dps) nsend[p]++; else nrecv[dps[0]]++; } Table send_data(nsend); Table recv_data(nrecv); /** Fill send_data **/ nsend = 0; for (auto i : newvr) { auto dps = GetDistantPNums(i); if (dps.Size() && rank < dps[0]) for(auto p : dps) send_data[p][nsend[p]++] = glob_vert[i]; } Array requests; for (int i = 0; i < comm.Size(); i++) { if (nsend[i]) requests.Append (comm.ISend (send_data[i], i, 200)); if (nrecv[i]) requests.Append (comm.IRecv (recv_data[i], i, 200)); } MyMPI_WaitAll (requests); Array cnt(comm.Size()); cnt = 0; for (auto i : newvr) { auto dps = GetDistantPNums(i); if (dps.Size() > 0 && dps[0] < comm.Rank()) { int master = comm.Size(); for (int j = 0; j < dps.Size(); j++) master = min (master, dps[j]); if (master != dps[0]) cout << "master not the first one !" << endl; glob_vert[i] = recv_data[master][cnt[master]++]; } } /* if (PointIndex::BASE==1) for (auto & i : glob_vert) i++; */ /* cout << "check ordering: " << endl; for (int i = 0; i < glob_vert.Size()-1; i++) if (glob_vert[i] > glob_vert[i+1]) cout << "wrong ordering" << endl; */ // reorder following global ordering: Array index0(glob_vert.Size()); for (int pi : Range(index0)) index0[pi] = pi; QuickSortI (FlatArray (glob_vert), index0); comm.Barrier(); for (int i = 0; i+1 < glob_vert.Size(); i++) if (glob_vert[index0[i]] > glob_vert[index0[i+1]]) cout << "wrong ordering" << endl; comm.Barrier(); if (rank != 0) { Array index(index0.Size()); for (int i = 0; i < index0.Size(); i++) index[i+PointIndex::BASE] = index0[i]+PointIndex::BASE; Array inv_index(index0.Size()); for (int i = 0; i < index0.Size(); i++) inv_index[index0[i]+PointIndex::BASE] = i+PointIndex::BASE; for (auto & el : mesh.VolumeElements()) for (PointIndex & pi : el.PNums()) pi = inv_index[pi]; for (auto & el : mesh.SurfaceElements()) for (PointIndex & pi : el.PNums()) pi = inv_index[pi]; for (auto & el : mesh.LineSegments()) for (PointIndex & pi : el.PNums()) pi = inv_index[pi]; // auto hpoints (mesh.Points()); Array hpoints { mesh.Points() }; for (PointIndex pi : Range(mesh.Points())) mesh.Points()[inv_index[pi]] = hpoints[pi]; if (mesh.mlbetweennodes.Size() == mesh.Points().Size()) { cout << "take care of multigrid table" << endl; NgArray,PointIndex::BASE> hml { mesh.mlbetweennodes }; for (PointIndex pi : Range(mesh.Points())) mesh.mlbetweennodes[inv_index[pi]] = hml[pi]; } // *testout << "index0 = " << endl << index0 << endl; // *testout << "loc2distvertold = " << endl; // for (auto i : Range(index0)) // *testout << "l " << i << " globi "<< glob_vert[i] << " dist = " << loc2distvert[i] << endl; DynamicTable oldtable(loc2distvert.Size()); for (size_t i = 0; i < loc2distvert.Size(); i++) for (auto val : loc2distvert[i]) oldtable.Add (i, val); loc2distvert = DynamicTable (oldtable.Size()); for (size_t i = 0; i < oldtable.Size(); i++) for (auto val : oldtable[index0[i]]) loc2distvert.Add (i, val); Array hglob_vert(glob_vert); for (int i = 0; i < index0.Size(); i++) glob_vert[i] = hglob_vert[index0[i]]; // *testout << "loc2distvertnew = " << endl; // for (auto i : Range(index0)) // *testout << "l " << i << " globi "<< glob_vert[i] << " dist = " << loc2distvert[i] << endl; } for (int i = 0; i+1 < glob_vert.Size(); i++) if (glob_vert[i] > glob_vert[i+1]) cout << "wrong ordering of globvert" << endl; // *testout << "new glob_vert = " << glob_vert << endl; } void ParallelMeshTopology :: SetDistantFaceNum (int dest, int locnum) { for ( int i = 0; i < loc2distface[locnum-1].Size(); i+=1 ) if ( loc2distface[locnum-1][i] == dest ) return; loc2distface.Add(locnum-1, dest); } void ParallelMeshTopology :: SetDistantPNum (int dest, int locnum) { for ( int i = 0; i < loc2distvert[locnum-1].Size(); i+=1 ) if ( loc2distvert[locnum-1][i] == dest ) return; loc2distvert.Add (locnum-1, dest); } void ParallelMeshTopology :: SetDistantEdgeNum (int dest, int locnum) { for ( int i = 0; i < loc2distedge[locnum-1].Size(); i+=1 ) if ( loc2distedge[locnum-1][i] == dest ) return; loc2distedge.Add (locnum-1, dest); } void ParallelMeshTopology :: SetNV_Loc2Glob (int anv) { glob_vert.SetSize(anv); glob_vert = -1; } void ParallelMeshTopology :: SetNV (int anv) { // glob_vert.SetSize(anv); // glob_vert = -1; // loc2distvert.ChangeSize (anv); DynamicTable oldtable(loc2distvert.Size()); for (size_t i = 0; i < loc2distvert.Size(); i++) for (auto val : loc2distvert[i]) oldtable.Add (i, val); loc2distvert = DynamicTable (anv); for (size_t i = 0; i < min(anv, oldtable.Size()); i++) for (auto val : oldtable[i]) loc2distvert.Add (i, val); } void ParallelMeshTopology :: SetNE ( int ane ) { glob_el.SetSize (ane); glob_el = -1; } void ParallelMeshTopology :: SetNSE ( int anse ) { glob_surfel.SetSize(anse); glob_surfel = -1; } void ParallelMeshTopology :: SetNSegm ( int anseg ) { glob_segm.SetSize (anseg); glob_segm = -1; } void ParallelMeshTopology :: UpdateCoarseGridGlobal () { // cout << "updatecoarsegridglobal called" << endl; if (id == 0) PrintMessage ( 3, "UPDATE GLOBAL COARSEGRID STARTS" ); int timer = NgProfiler::CreateTimer ("UpdateCoarseGridGlobal"); NgProfiler::RegionTimer reg(timer); *testout << "ParallelMeshTopology :: UpdateCoarseGridGlobal" << endl; const MeshTopology & topology = mesh.GetTopology(); auto comm = mesh.GetCommunicator(); if ( id == 0 ) { NgArray*> sendarrays(ntasks); for (int dest = 1; dest < ntasks; dest++) sendarrays[dest] = new NgArray; NgArray edges, faces; for (int el = 1; el <= mesh.GetNE(); el++) { topology.GetElementFaces (el, faces); topology.GetElementEdges (el, edges); const Element & volel = mesh.VolumeElement (el); // NgArray & sendarray = *sendarrays[volel.GetPartition()]; NgArray & sendarray = *sendarrays[mesh.vol_partition[el-1]]; for ( int i = 0; i < edges.Size(); i++ ) sendarray.Append (edges[i]); for ( int i = 0; i < faces.Size(); i++ ) sendarray.Append (faces[i]); } for (int el = 1; el <= mesh.GetNSE(); el++) { topology.GetSurfaceElementEdges (el, edges); const Element2d & surfel = mesh.SurfaceElement (el); // NgArray & sendarray = *sendarrays[surfel.GetPartition()]; NgArray & sendarray = *sendarrays[mesh.surf_partition[el-1]]; for ( int i = 0; i < edges.Size(); i++ ) sendarray.Append (edges[i]); sendarray.Append (topology.GetSurfaceElementFace (el)); } Array sendrequests; for (int dest = 1; dest < ntasks; dest++) // sendrequests.Append (MyMPI_ISend (*sendarrays[dest], dest, MPI_TAG_MESH+10, comm)); sendrequests.Append (comm.ISend (FlatArray(*sendarrays[dest]), dest, MPI_TAG_MESH+10)); MyMPI_WaitAll (sendrequests); for (int dest = 1; dest < ntasks; dest++) delete sendarrays[dest]; } else { NgArray recvarray; MyMPI_Recv (recvarray, 0, MPI_TAG_MESH+10, comm); int ii = 0; NgArray faces, edges; for (int volel = 1; volel <= mesh.GetNE(); volel++) { topology.GetElementEdges ( volel, edges); for ( int i = 0; i < edges.Size(); i++) SetLoc2Glob_Edge ( edges[i], recvarray[ii++]); topology.GetElementFaces( volel, faces); for ( int i = 0; i < faces.Size(); i++) SetLoc2Glob_Face ( faces[i], recvarray[ii++]); } for (int surfel = 1; surfel <= mesh.GetNSE(); surfel++) { topology.GetSurfaceElementEdges (surfel, edges); for (int i = 0; i < edges.Size(); i++) SetLoc2Glob_Edge (edges[i], recvarray[ii++]); int face = topology.GetSurfaceElementFace (surfel); SetLoc2Glob_Face ( face, recvarray[ii++]); } } is_updated = true; } void ParallelMeshTopology :: IdentifyVerticesAfterRefinement() { static Timer t("ParallelTopology::UpdateCoarseGrid"); RegionTimer r(t); // cout << "UpdateCoarseGrid" << endl; // if (is_updated) return; NgMPI_Comm comm = mesh.GetCommunicator(); int id = comm.Rank(); int ntasks = comm.Size(); if (ntasks == 1) return; Reset(); static int timer = NgProfiler::CreateTimer ("UpdateCoarseGrid"); NgProfiler::RegionTimer reg(timer); (*testout) << "UPDATE COARSE GRID PARALLEL TOPOLOGY " << endl; if (id == 0) PrintMessage (1, "update parallel topology"); // UpdateCoarseGridGlobal(); // MPI_Barrier (MPI_COMM_WORLD); MPI_Group MPI_GROUP_comm; MPI_Group MPI_LocalGroup; MPI_Comm MPI_LocalComm1; int process_ranks[] = { 0 }; MPI_Comm_group (comm, &MPI_GROUP_comm); MPI_Group_excl (MPI_GROUP_comm, 1, process_ranks, &MPI_LocalGroup); MPI_Comm_create (comm, MPI_LocalGroup, &MPI_LocalComm1); if (id == 0) { // SetNV(0); // EnumeratePointsGlobally(); return; } NgMPI_Comm MPI_LocalComm(MPI_LocalComm1); const MeshTopology & topology = mesh.GetTopology(); NgArray cnt_send(ntasks-1); // update new vertices after mesh-refinement if (mesh.mlbetweennodes.Size() > 0) { // *testout << "have to identify new vertices, nv = " << mesh.GetNV() << endl; // cout << "UpdateCoarseGrid - vertices" << endl; int newnv = mesh.mlbetweennodes.Size(); // loc2distvert.ChangeSize(mesh.mlbetweennodes.Size()); DynamicTable oldtable(loc2distvert.Size()); for (size_t i = 0; i < loc2distvert.Size(); i++) for (auto val : loc2distvert[i]) oldtable.Add (i, val); loc2distvert = DynamicTable (mesh.mlbetweennodes.Size()); for (size_t i = 0; i < min(loc2distvert.Size(), oldtable.Size()); i++) for (auto val : oldtable[i]) loc2distvert.Add (i, val); *testout << "extended loc2distver = " << endl << loc2distvert << endl; /* for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) for (int dest = 1; dest < ntasks; dest++) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) SetDistantPNum(dest, pi); } */ bool changed = true; while (changed) { changed = false; // build exchange vertices cnt_send = 0; for (PointIndex pi : mesh.Points().Range()) for (int dist : GetDistantPNums(pi-PointIndex::BASE)) cnt_send[dist-1]++; TABLE dest2vert(cnt_send); for (PointIndex pi : mesh.Points().Range()) for (int dist : GetDistantPNums(pi-PointIndex::BASE)) dest2vert.Add (dist-1, pi); for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) // for (int dest = 1; dest < ntasks; dest++) for (int dest : GetDistantPNums(v1-PointIndex::BASE)) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) cnt_send[dest-1]++; } TABLE dest2pair(cnt_send); // for (int dest = 1; dest < ntasks; dest++) for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) for (int dest : GetDistantPNums(v1-PointIndex::BASE)) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) dest2pair.Add (dest-1, pi); } cnt_send = 0; int v1, v2; for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) for (int dest : GetDistantPNums(v1-PointIndex::BASE)) if (IsExchangeVert(dest, v2)) cnt_send[dest-1]+=2; } TABLE send_verts(cnt_send); NgArray loc2exchange(mesh.GetNV()); for (int dest = 1; dest < ntasks; dest++) if (dest != id) { loc2exchange = -1; int cnt = 0; /* for (PointIndex pi : mesh.Points().Range()) if (IsExchangeVert(dest, pi)) loc2exchange[pi] = cnt++; */ for (PointIndex pi : dest2vert[dest-1]) loc2exchange[pi] = cnt++; // for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) for (PointIndex pi : dest2pair[dest-1]) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) { send_verts.Add (dest-1, loc2exchange[v1]); send_verts.Add (dest-1, loc2exchange[v2]); } } } TABLE recv_verts(ntasks-1); MyMPI_ExchangeTable (send_verts, recv_verts, MPI_TAG_MESH+9, MPI_LocalComm); for (int dest = 1; dest < ntasks; dest++) if (dest != id) { loc2exchange = -1; int cnt = 0; /* for (PointIndex pi : mesh.Points().Range()) if (IsExchangeVert(dest, pi)) loc2exchange[pi] = cnt++; */ for (PointIndex pi : dest2vert[dest-1]) loc2exchange[pi] = cnt++; NgFlatArray recvarray = recv_verts[dest-1]; for (int ii = 0; ii < recvarray.Size(); ii+=2) for (PointIndex pi : dest2pair[dest-1]) // for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) { INDEX_2 re(recvarray[ii], recvarray[ii+1]); INDEX_2 es(loc2exchange[v1], loc2exchange[v2]); if (es == re && !IsExchangeVert(dest, pi)) { SetDistantPNum(dest, pi); changed = true; } } } } } } NgArray sendarray, recvarray; // cout << "UpdateCoarseGrid - edges" << endl; // static int timerv = NgProfiler::CreateTimer ("UpdateCoarseGrid - ex vertices"); static int timere = NgProfiler::CreateTimer ("UpdateCoarseGrid - ex edges"); static int timerf = NgProfiler::CreateTimer ("UpdateCoarseGrid - ex faces"); NgProfiler::StartTimer (timere); // build exchange vertices cnt_send = 0; for (PointIndex pi : mesh.Points().Range()) for (int dist : GetDistantPNums(pi-PointIndex::BASE)) cnt_send[dist-1]++; TABLE dest2vert(cnt_send); for (PointIndex pi : mesh.Points().Range()) for (int dist : GetDistantPNums(pi-PointIndex::BASE)) dest2vert.Add (dist-1, pi); MPI_Group_free(&MPI_LocalGroup); // MPI_Comm_free(&MPI_LocalComm); } void ParallelMeshTopology :: UpdateCoarseGrid () { static Timer t("ParallelTopology::UpdateCoarseGrid"); RegionTimer r(t); // cout << "UpdateCoarseGrid" << endl; // if (is_updated) return; NgMPI_Comm comm = mesh.GetCommunicator(); int id = comm.Rank(); int ntasks = comm.Size(); if (ntasks == 1) return; Reset(); static int timer = NgProfiler::CreateTimer ("UpdateCoarseGrid"); NgProfiler::RegionTimer reg(timer); (*testout) << "UPDATE COARSE GRID PARALLEL TOPOLOGY " << endl; if (id == 0) PrintMessage (1, "update parallel topology"); // UpdateCoarseGridGlobal(); // MPI_Barrier (MPI_COMM_WORLD); MPI_Group MPI_GROUP_comm; MPI_Group MPI_LocalGroup; MPI_Comm MPI_LocalComm; int process_ranks[] = { 0 }; MPI_Comm_group (comm, &MPI_GROUP_comm); MPI_Group_excl (MPI_GROUP_comm, 1, process_ranks, &MPI_LocalGroup); MPI_Comm_create (comm, MPI_LocalGroup, &MPI_LocalComm); if (id == 0) { // SetNV(0); // EnumeratePointsGlobally(); return; } const MeshTopology & topology = mesh.GetTopology(); NgArray cnt_send(ntasks-1); #ifdef NONE // update new vertices after mesh-refinement if (mesh.mlbetweennodes.Size() > 0) { // cout << "UpdateCoarseGrid - vertices" << endl; int newnv = mesh.mlbetweennodes.Size(); // loc2distvert.ChangeSize(mesh.mlbetweennodes.Size()); DynamicTable oldtable(loc2distvert.Size()); for (size_t i = 0; i < loc2distvert.Size(); i++) for (auto val : loc2distvert[i]) oldtable.Add (i, val); loc2distvert = DynamicTable (mesh.mlbetweennodes.Size()); for (size_t i = 0; i < min(loc2distvert.Size(), oldtable.Size()); i++) for (auto val : oldtable[i]) loc2distvert.Add (i, val); /* for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) for (int dest = 1; dest < ntasks; dest++) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) SetDistantPNum(dest, pi); } */ bool changed = true; while (changed) { changed = false; // build exchange vertices cnt_send = 0; for (PointIndex pi : mesh.Points().Range()) for (int dist : GetDistantPNums(pi-PointIndex::BASE)) cnt_send[dist-1]++; TABLE dest2vert(cnt_send); for (PointIndex pi : mesh.Points().Range()) for (int dist : GetDistantPNums(pi-PointIndex::BASE)) dest2vert.Add (dist-1, pi); for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) // for (int dest = 1; dest < ntasks; dest++) for (int dest : GetDistantPNums(v1-PointIndex::BASE)) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) cnt_send[dest-1]++; } TABLE dest2pair(cnt_send); // for (int dest = 1; dest < ntasks; dest++) for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) for (int dest : GetDistantPNums(v1-PointIndex::BASE)) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) dest2pair.Add (dest-1, pi); } cnt_send = 0; int v1, v2; for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) for (int dest : GetDistantPNums(v1-PointIndex::BASE)) if (IsExchangeVert(dest, v2)) cnt_send[dest-1]+=2; } TABLE send_verts(cnt_send); NgArray loc2exchange(mesh.GetNV()); for (int dest = 1; dest < ntasks; dest++) if (dest != id) { loc2exchange = -1; int cnt = 0; /* for (PointIndex pi : mesh.Points().Range()) if (IsExchangeVert(dest, pi)) loc2exchange[pi] = cnt++; */ for (PointIndex pi : dest2vert[dest-1]) loc2exchange[pi] = cnt++; // for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) for (PointIndex pi : dest2pair[dest-1]) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) { send_verts.Add (dest-1, loc2exchange[v1]); send_verts.Add (dest-1, loc2exchange[v2]); } } } TABLE recv_verts(ntasks-1); MyMPI_ExchangeTable (send_verts, recv_verts, MPI_TAG_MESH+9, MPI_LocalComm); for (int dest = 1; dest < ntasks; dest++) if (dest != id) { loc2exchange = -1; int cnt = 0; /* for (PointIndex pi : mesh.Points().Range()) if (IsExchangeVert(dest, pi)) loc2exchange[pi] = cnt++; */ for (PointIndex pi : dest2vert[dest-1]) loc2exchange[pi] = cnt++; NgFlatArray recvarray = recv_verts[dest-1]; for (int ii = 0; ii < recvarray.Size(); ii+=2) for (PointIndex pi : dest2pair[dest-1]) // for (PointIndex pi = PointIndex::BASE; pi < newnv+PointIndex::BASE; pi++) { PointIndex v1 = mesh.mlbetweennodes[pi][0]; PointIndex v2 = mesh.mlbetweennodes[pi][1]; if (mesh.mlbetweennodes[pi][0] != PointIndex::BASE-1) { INDEX_2 re(recvarray[ii], recvarray[ii+1]); INDEX_2 es(loc2exchange[v1], loc2exchange[v2]); if (es == re && !IsExchangeVert(dest, pi)) { SetDistantPNum(dest, pi); changed = true; } } } } } } #endif NgArray sendarray, recvarray; // cout << "UpdateCoarseGrid - edges" << endl; // static int timerv = NgProfiler::CreateTimer ("UpdateCoarseGrid - ex vertices"); static int timere = NgProfiler::CreateTimer ("UpdateCoarseGrid - ex edges"); static int timerf = NgProfiler::CreateTimer ("UpdateCoarseGrid - ex faces"); NgProfiler::StartTimer (timere); int nfa = topology . GetNFaces(); int ned = topology . GetNEdges(); // build exchange vertices cnt_send = 0; for (PointIndex pi : mesh.Points().Range()) for (int dist : GetDistantPNums(pi-PointIndex::BASE)) cnt_send[dist-1]++; TABLE dest2vert(cnt_send); for (PointIndex pi : mesh.Points().Range()) for (int dist : GetDistantPNums(pi-PointIndex::BASE)) dest2vert.Add (dist-1, pi); // exchange edges cnt_send = 0; int v1, v2; for (int edge = 1; edge <= ned; edge++) { topology.GetEdgeVertices (edge, v1, v2); for (int dest = 1; dest < ntasks; dest++) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) cnt_send[dest-1]+=1; } TABLE dest2edge(cnt_send); for (int & v : cnt_send) v *= 2; TABLE send_edges(cnt_send); for (int edge = 1; edge <= ned; edge++) { topology.GetEdgeVertices (edge, v1, v2); for (int dest = 1; dest < ntasks; dest++) if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) dest2edge.Add (dest-1, edge); } NgArray loc2exchange(mesh.GetNV()); for (int dest = 1; dest < ntasks; dest++) { loc2exchange = -1; int cnt = 0; for (PointIndex pi : dest2vert[dest-1]) loc2exchange[pi] = cnt++; for (int edge : dest2edge[dest-1]) { topology.GetEdgeVertices (edge, v1, v2); if (IsExchangeVert (dest, v1) && IsExchangeVert (dest, v2)) { send_edges.Add (dest-1, loc2exchange[v1]); send_edges.Add (dest-1, loc2exchange[v2]); } } } // cout << "UpdateCoarseGrid - edges mpi-exchange" << endl; TABLE recv_edges(ntasks-1); MyMPI_ExchangeTable (send_edges, recv_edges, MPI_TAG_MESH+9, MPI_LocalComm); // cout << "UpdateCoarseGrid - edges mpi-exchange done" << endl; for (int dest = 1; dest < ntasks; dest++) { auto ex2loc = dest2vert[dest-1]; if (ex2loc.Size() == 0) continue; INDEX_2_CLOSED_HASHTABLE vert2edge(2*dest2edge[dest-1].Size()+10); for (int edge : dest2edge[dest-1]) { topology.GetEdgeVertices (edge, v1, v2); vert2edge.Set(INDEX_2(v1,v2), edge); } NgFlatArray recvarray = recv_edges[dest-1]; for (int ii = 0; ii < recvarray.Size(); ii+=2) { INDEX_2 re(ex2loc[recvarray[ii]], ex2loc[recvarray[ii+1]]); if (vert2edge.Used(re)) SetDistantEdgeNum(dest, vert2edge.Get(re)); } } NgProfiler::StopTimer (timere); // cout << "UpdateCoarseGrid - faces" << endl; if (mesh.GetDimension() == 3) { NgProfiler::StartTimer (timerf); NgArray verts; // exchange faces cnt_send = 0; for (int face = 1; face <= nfa; face++) { topology.GetFaceVertices (face, verts); for (int dest = 1; dest < ntasks; dest++) if (dest != id) if (IsExchangeVert (dest, verts[0]) && IsExchangeVert (dest, verts[1]) && IsExchangeVert (dest, verts[2])) cnt_send[dest-1]++; } TABLE dest2face(cnt_send); for (int face = 1; face <= nfa; face++) { topology.GetFaceVertices (face, verts); for (int dest = 1; dest < ntasks; dest++) if (dest != id) if (IsExchangeVert (dest, verts[0]) && IsExchangeVert (dest, verts[1]) && IsExchangeVert (dest, verts[2])) dest2face.Add(dest-1, face); } for (int & c : cnt_send) c*=3; TABLE send_faces(cnt_send); NgArray loc2exchange(mesh.GetNV()); for (int dest = 1; dest < ntasks; dest++) if (dest != id) { if (dest2vert[dest-1].Size() == 0) continue; loc2exchange = -1; int cnt = 0; for (PointIndex pi : dest2vert[dest-1]) loc2exchange[pi] = cnt++; for (int face : dest2face[dest-1]) { topology.GetFaceVertices (face, verts); if (IsExchangeVert (dest, verts[0]) && IsExchangeVert (dest, verts[1]) && IsExchangeVert (dest, verts[2])) { send_faces.Add (dest-1, loc2exchange[verts[0]]); send_faces.Add (dest-1, loc2exchange[verts[1]]); send_faces.Add (dest-1, loc2exchange[verts[2]]); } } } // cout << "UpdateCoarseGrid - faces mpi-exchange" << endl; TABLE recv_faces(ntasks-1); MyMPI_ExchangeTable (send_faces, recv_faces, MPI_TAG_MESH+9, MPI_LocalComm); // cout << "UpdateCoarseGrid - faces mpi-exchange done" << endl; for (int dest = 1; dest < ntasks; dest++) { auto ex2loc = dest2vert[dest-1]; if (ex2loc.Size() == 0) continue; INDEX_3_CLOSED_HASHTABLE vert2face(2*dest2face[dest-1].Size()+10); for (int face : dest2face[dest-1]) { topology.GetFaceVertices (face, verts); vert2face.Set(INDEX_3(verts[0], verts[1], verts[2]), face); } NgFlatArray recvarray = recv_faces[dest-1]; for (int ii = 0; ii < recvarray.Size(); ii+=3) { INDEX_3 re(ex2loc[recvarray[ii]], ex2loc[recvarray[ii+1]], ex2loc[recvarray[ii+2]]); if (vert2face.Used(re)) SetDistantFaceNum(dest, vert2face.Get(re)); } } NgProfiler::StopTimer (timerf); } // cout << "UpdateCoarseGrid - done" << endl; // EnumeratePointsGlobally(); is_updated = true; MPI_Group_free(&MPI_LocalGroup); MPI_Comm_free(&MPI_LocalComm); } } #endif