#include "analysis/analysis.h" #include "ProfilerApp.h" #include template inline TYPE* getPtr( std::vector& x ) { return x.empty() ? NULL:&x[0]; } template inline const TYPE* getPtr( const std::vector& x ) { return x.empty() ? NULL:&x[0]; } /****************************************************************** * Compute the blobs * ******************************************************************/ int ComputeBlob( const Array& isPhase, BlobIDArray& LocalBlobID, bool periodic, int start_id ) { PROFILE_START("ComputeBlob",1); ASSERT(isPhase.size()==LocalBlobID.size()); const int Nx = isPhase.size(0); // maxima for the meshes const int Ny = isPhase.size(1); const int Nz = isPhase.size(2); std::vector map; map.reserve(128); // Get the list of neighbors we need to check int N_neighbors = 0; int d[26][3]; bool include_corners = false; // Do we need to include cells that only touch at a corder/edge if ( include_corners ) { // Include corners/edges as neighbors, check all cells N_neighbors = 26; const int tmp[26][3] = {{1,0,0},{-1,0,0},{0,1,0},{0,-1,0},{0,0,1},{0,0,-1}, {1,1,0},{1,-1,0},{-1,1,0},{-1,-1,0},{1,0,1},{-1,0,1}, {1,0,-1},{-1,0,-1},{0,1,1},{0,-1,1},{0,1,-1},{0,-1,-1}, {1,1,1},{1,1,-1},{1,-1,1},{1,-1,-1},{-1,1,1},{-1,1,-1}, {-1,-1,1},{-1,-1,-1}}; // directions to neighbors memcpy(d,tmp,sizeof(tmp)); } else { // Do not include corners/edges as neighbors if ( periodic ) { // Include all neighbors for periodic problems N_neighbors = 6; const int tmp[6][3] = {{1,0,0},{-1,0,0},{0,1,0},{0,-1,0},{0,0,1},{0,0,-1}}; // directions to neighbors memcpy(d,tmp,sizeof(tmp)); } else { // We only need to include the lower neighbors for non-periodic problems N_neighbors = 3; const int tmp[3][3] = {{-1,0,0},{0,-1,0},{0,0,-1}}; // directions to neighbors memcpy(d,tmp,sizeof(tmp)); } } // Loop through all the points int last = start_id-1; std::vector neighbor_ids; neighbor_ids.reserve(N_neighbors); const bool *isPhasePtr = isPhase.get(); BlobIDType *LocalBlobIDPtr = LocalBlobID.get(); for (int z=0; zNx-1 ? 0:x2; y2 = y2<0 ? Ny-1:y2; // Periodic BC for x y2 = y2>Ny-1 ? 0:y2; z2 = z2<0 ? Nz-1:z2; // Periodic BC for x z2 = z2>Nz-1 ? 0:z2; } else { if ( x2<0 || x2>=Nx || y2<0 || y2>=Ny || z2<0 || z2>=Nz ) continue; } // Check if a neighbor has a known blob id size_t index2 = x2 + y2*Nx + z2*Nx*Ny; int id = LocalBlobIDPtr[index2]; if ( !isPhasePtr[index2] || id<0 ) continue; neighbor_ids[N_list] = id; N_list++; } if ( N_list==0 ) { // No neighbors with a blob id, create a new one LocalBlobIDPtr[index] = last+1; map.push_back(last+1); last++; } else if ( N_list==1 ) { // We have one neighbor LocalBlobIDPtr[index] = neighbor_ids[0]; } else { // We have multiple neighbors int id = neighbor_ids[0]; for (int i=1; i isPhase(Nx,Ny,Nz); memset(isPhase.get(),0,Nx*Ny*Nz*sizeof(bool)); for (size_t i=0; ivF && SignDist(i)>vS ) isPhase(i) = true; } } int nblobs = ComputeBlob( isPhase, LocalBlobID, periodic, 0 ); PROFILE_STOP("ComputeLocalBlobIDs"); return nblobs; } int ComputeLocalPhaseComponent(const IntArray &PhaseID, int VALUE, BlobIDArray &ComponentLabel, bool periodic ) { PROFILE_START("ComputeLocalPhaseComponent"); size_t Nx = PhaseID.size(0); size_t Ny = PhaseID.size(1); size_t Nz = PhaseID.size(2); size_t N = Nx*Ny*Nz; // Compute the local blob ids ComponentLabel.resize(Nx,Ny,Nz); Array isPhase(Nx,Ny,Nz); for (size_t i=0; i= 0 ) local_size[id] += 1; max_id = std::max(max_id,id); } } } ASSERT(max_id > map1(N_blobs); int N_blobs2 = 0; for (int i=0; i 0 ) N_blobs2++; } std::sort( map1.begin(), map1.end() ); std::vector map2(N_blobs,-1); for (int i=0; i= 0 ) ID(i) = map2[ID(i)]; } delete [] local_size; delete [] global_size; PROFILE_STOP("ReorderBlobIDs2",1); return N_blobs2; } void ReorderBlobIDs( BlobIDArray& ID ) { PROFILE_START("ReorderBlobIDs"); int tmp = ID.max()+1; int N_blobs = 0; MPI_Allreduce(&tmp,&N_blobs,1,MPI_INT,MPI_MAX,MPI_COMM_WORLD); ReorderBlobIDs2(ID,N_blobs,1,1,1); PROFILE_STOP("ReorderBlobIDs"); } /****************************************************************** * Compute the global blob ids * ******************************************************************/ struct global_id_info_struct { int64_t new_id; std::set remote_ids; }; // Send the local ids and their new value to all neighbors static void updateRemoteIds( const std::map& map, const std::vector& neighbors, int N_send, const std::vector& N_recv, int64_t *send_buf, std::vector& recv_buf, std::map& remote_map ) { std::vector send_req(neighbors.size()); std::vector recv_req(neighbors.size()); std::vector status(neighbors.size()); std::map::const_iterator it = map.begin(); ASSERT(N_send==(int)map.size()); for (size_t i=0; ifirst; send_buf[2*i+1] = it->second.new_id; } for (size_t i=0; ifirst] = it->second.new_id; } for (size_t i=0; i& remote_map, std::map& map ) { bool changed = false; std::map::iterator it; for (it=map.begin(); it!=map.end(); ++it) { int64_t id = it->second.new_id; std::set::const_iterator it2; for (it2=it->second.remote_ids.begin(); it2!=it->second.remote_ids.end(); ++it2) { int64_t id2 = remote_map.find(*it2)->second; id = std::min(id,id2); } changed = changed || it->second.new_id!=id; it->second.new_id = id; } return changed; } static int LocalToGlobalIDs( int nx, int ny, int nz, const RankInfoStruct& rank_info, int nblobs, BlobIDArray& IDs ) { PROFILE_START("LocalToGlobalIDs",1); const int rank = rank_info.rank[1][1][1]; int nprocs; MPI_Comm_size(MPI_COMM_WORLD,&nprocs); const int ngx = (IDs.size(0)-nx)/2; const int ngy = (IDs.size(1)-ny)/2; const int ngz = (IDs.size(2)-nz)/2; // Get the number of blobs for each rank std::vector N_blobs(nprocs,0); PROFILE_START("LocalToGlobalIDs-Allgather",1); MPI_Allgather(&nblobs,1,MPI_INT,getPtr(N_blobs),1,MPI_INT,MPI_COMM_WORLD); PROFILE_STOP("LocalToGlobalIDs-Allgather",1); int64_t N_blobs_tot = 0; int offset = 0; for (int i=0; i= 0 ) IDs(i) += offset; } const BlobIDArray LocalIDs = IDs; // Copy the ids and get the neighbors through the halos fillHalo fillData(rank_info,nx,ny,nz,1,1,1,0,1,true,true,true); fillData.fill(IDs); // Create a list of all neighbor ranks (excluding self) std::vector neighbors; neighbors.push_back( rank_info.rank[0][1][1] ); neighbors.push_back( rank_info.rank[2][1][1] ); neighbors.push_back( rank_info.rank[1][0][1] ); neighbors.push_back( rank_info.rank[1][2][1] ); neighbors.push_back( rank_info.rank[1][1][0] ); neighbors.push_back( rank_info.rank[1][1][2] ); std::unique(neighbors.begin(),neighbors.end()); if ( std::find(neighbors.begin(),neighbors.end(),rank) != neighbors.end() ) neighbors.erase(std::find(neighbors.begin(),neighbors.end(),rank)); // Create a map of all local ids to the neighbor ids std::map map; std::set local; for (size_t i=0; i=0 ) { local.insert(LocalIDs(i)); if ( LocalIDs(i)!=IDs(i) && IDs(i)>=0 ) map[LocalIDs(i)].remote_ids.insert(IDs(i)); } } std::map::iterator it; for (it=map.begin(); it!=map.end(); ++it) { it->second.new_id = it->first; local.erase(it->first); } // Get the number of ids we will recieve from each rank int N_send = map.size(); std::vector N_recv(neighbors.size(),0); std::vector send_req(neighbors.size()); std::vector recv_req(neighbors.size()); std::vector status(neighbors.size()); for (size_t i=0; i recv_buf(neighbors.size()); for (size_t i=0; i remote_map; for (it=map.begin(); it!=map.end(); ++it) { int64_t id = it->first; std::set::const_iterator it2; for (it2=it->second.remote_ids.begin(); it2!=it->second.remote_ids.end(); ++it2) { int64_t id2 = *it2; id = std::min(id,id2); remote_map.insert(std::pair(id2,id2)); } it->second.new_id = id; } // Iterate until we are done int iteration = 1; PROFILE_START("LocalToGlobalIDs-loop",1); while ( 1 ) { iteration++; // Send the local ids and their new value to all neighbors updateRemoteIds( map, neighbors, N_send, N_recv,send_buf, recv_buf, remote_map ); // Compute a new local id for each local id bool changed = updateLocalIds( remote_map, map ); // Check if we are finished int test = changed ? 1:0; int result = 0; MPI_Allreduce(&test,&result,1,MPI_INT,MPI_SUM,MPI_COMM_WORLD); if ( result==0 ) break; } PROFILE_STOP("LocalToGlobalIDs-loop",1); // Relabel the ids std::vector final_map(nblobs,-1); for (it=map.begin(); it!=map.end(); ++it) final_map[it->first-offset] = it->second.new_id; for (std::set::const_iterator it2=local.begin(); it2!=local.end(); ++it2) final_map[*it2-offset] = *it2; for (size_t i=0; i=0); for (size_t k=ngz; k= 0 ) IDs(i,j,k) = final_map[id-offset]; } } } // Fill the ghosts fillHalo fillData2(rank_info,nx,ny,nz,1,1,1,0,1,true,true,true); fillData2.fill(IDs); // Reorder based on size (and compress the id space int N_blobs_global = ReorderBlobIDs2(IDs,N_blobs_tot,ngx,ngy,ngz); // Finished delete [] send_buf; for (size_t i=0; i inline MPI_Datatype getMPIType(); template<> inline MPI_Datatype getMPIType() { return MPI_INT; } template<> inline MPI_Datatype getMPIType() { if ( sizeof(int64_t)==sizeof(long int) ) return MPI_LONG; else if ( sizeof(int64_t)==sizeof(double) ) return MPI_DOUBLE; } template void gatherSet( std::set& set ) { int nprocs; MPI_Comm_size(MPI_COMM_WORLD,&nprocs); MPI_Datatype type = getMPIType(); std::vector send_data(set.begin(),set.end()); int send_count = send_data.size(); std::vector recv_count(nprocs,0), recv_disp(nprocs,0); MPI_Allgather(&send_count,1,MPI_INT,getPtr(recv_count),1,MPI_INT,MPI_COMM_WORLD); for (int i=1; i recv_data(recv_disp[nprocs-1]+recv_count[nprocs-1]); MPI_Allgatherv(getPtr(send_data),send_count,type, getPtr(recv_data),getPtr(recv_count),getPtr(recv_disp),type, MPI_COMM_WORLD); for (size_t i=0; i void gatherSrcIDMap( std::map >& src_map ) { int nprocs; MPI_Comm_size(MPI_COMM_WORLD,&nprocs); MPI_Datatype type = getMPIType(); std::vector send_data; typename std::map >::const_iterator it; for (it=src_map.begin(); it!=src_map.end(); ++it) { int id = it->first; const std::set& src_ids = it->second; send_data.push_back(id); send_data.push_back(src_ids.size()); typename std::set::const_iterator it2; for (it2=src_ids.begin(); it2!=src_ids.end(); ++it2) send_data.push_back(*it2); } int send_count = send_data.size(); std::vector recv_count(nprocs,0), recv_disp(nprocs,0); MPI_Allgather(&send_count,1,MPI_INT,getPtr(recv_count),1,MPI_INT,MPI_COMM_WORLD); for (int i=1; i recv_data(recv_disp[nprocs-1]+recv_count[nprocs-1]); MPI_Allgatherv(getPtr(send_data),send_count,type, getPtr(recv_data),getPtr(recv_count),getPtr(recv_disp),type, MPI_COMM_WORLD); size_t i=0; while ( i < recv_data.size() ) { int id = recv_data[i]; int count = recv_data[i+1]; i += 2; std::set& src_ids = src_map[id]; for (int j=0; j >& src_map, std::map >& dst_map, std::set& src, std::set& dst ) { src.insert(src_id); const std::set& dst_ids = dst_map[src_id]; for (std::set::const_iterator it=dst_ids.begin(); it!=dst_ids.end(); ++it) { if ( dst.find(*it)==dst.end() ) addSrcDstIDs(*it,dst_map,src_map,dst,src); } } ID_map_struct computeIDMap( const BlobIDArray& ID1, const BlobIDArray& ID2 ) { ASSERT(ID1.size()==ID2.size()); PROFILE_START("computeIDMap"); // Get a global list of all src/dst ids and the src map for each local blob std::set src_set, dst_set; std::map > src_map; // Map of the src ids for each dst id for (size_t i=0; i=0 ) src_set.insert(ID1(i)); if ( ID2(i)>=0 ) dst_set.insert(ID2(i)); if ( ID2(i)>=0 && ID1(i)>=0 ) { std::set& src_ids = src_map[ID2(i)]; src_ids.insert(ID1(i)); } } // Communicate the src/dst ids and src id map to all processors and reduce gatherSet( src_set ); gatherSet( dst_set ); gatherSrcIDMap( src_map ); // Compute the dst id map std::map > dst_map; // Map of the dst ids for each src id for (std::map >::const_iterator it=src_map.begin(); it!=src_map.end(); ++it) { BlobIDType id = it->first; const std::set& src_ids = it->second; for (std::set::const_iterator it2=src_ids.begin(); it2!=src_ids.end(); ++it2) { std::set& dst_ids = dst_map[*it2]; dst_ids.insert(id); } } // Perform the mapping of ids ID_map_struct id_map; // Find new blobs for (std::set::const_iterator it=dst_set.begin(); it!=dst_set.end(); ++it) { if ( src_map.find(*it)==src_map.end() ) id_map.created.push_back(*it); } // Fine blobs that disappeared for (std::set::const_iterator it=src_set.begin(); it!=src_set.end(); ++it) { if ( dst_map.find(*it)==dst_map.end() ) id_map.destroyed.push_back(*it); } // Find blobs with a 1-to-1 mapping std::vector dst_list; dst_list.reserve(src_map.size()); for (std::map >::const_iterator it=src_map.begin(); it!=src_map.end(); ++it) dst_list.push_back(it->first); for (size_t i=0; i& src_ids = src_map[dst_id]; if ( src_ids.size()==1 ) { int src_id = *src_ids.begin(); const std::set& dst_ids = dst_map[src_id]; if ( dst_ids.size()==1 ) { ASSERT(*dst_ids.begin()==dst_id); src_map.erase(dst_id); dst_map.erase(src_id); id_map.src_dst.push_back(std::pair(src_id,dst_id)); } } } // Handle merge/splits while ( !dst_map.empty() ) { // Get a lit of the src-dst ids std::set src, dst; addSrcDstIDs( dst_map.begin()->first, src_map, dst_map, src, dst ); for (std::set::const_iterator it=src.begin(); it!=src.end(); ++it) dst_map.erase(*it); for (std::set::const_iterator it=dst.begin(); it!=dst.end(); ++it) src_map.erase(*it); if ( src.size()==1 ) { // Bubble split id_map.split.push_back( BlobIDSplitStruct(*src.begin(),std::vector(dst.begin(),dst.end())) ); } else if ( dst.size()==1 ) { // Bubble merge id_map.merge.push_back( BlobIDMergeStruct(std::vector(src.begin(),src.end()),*dst.begin()) ); } else { // Bubble split/merge id_map.merge_split.push_back( BlobIDMergeSplitStruct( std::vector(src.begin(),src.end()), std::vector(dst.begin(),dst.end() ) ) ); } } ASSERT(src_map.empty()); PROFILE_STOP("computeIDMap"); return id_map; } /****************************************************************** * Renumber the ids * ******************************************************************/ void getNewIDs( ID_map_struct& map, BlobIDType& id_max, std::vector& new_ids ) { new_ids.resize(0); // Renumber the ids that map directly for (size_t i=0; i(id1+1) ) new_ids.resize(id1+1,-1); new_ids[id1] = id2; } // Renumber the created blobs to create new ids for (size_t i=0; i(id1+1) ) new_ids.resize(id1+1,-1); new_ids[id1] = id2; } // Renumber the blob splits to create new ids for (size_t i=0; i(id1+1) ) new_ids.resize(id1+1,-1); new_ids[id1] = id2; } } // Renumber the blob merges to create a new id for (size_t i=0; i(id1+1) ) new_ids.resize(id1+1,-1); new_ids[id1] = id2; } // Renumber the blob merge/splits to create new ids for (size_t i=0; i(id1+1) ) new_ids.resize(id1+1,-1); new_ids[id1] = id2; } } } void renumberIDs( const std::vector& new_ids, BlobIDArray& IDs ) { size_t N = IDs.length(); BlobIDType* ids = IDs.get(); for (size_t i=0; i=0 ) ids[i] = new_ids[id]; } } /****************************************************************** * Write the id map for the given timestep * ******************************************************************/ void writeIDMap( const ID_map_struct& map, long long int timestep, const std::string& filename ) { int rank; MPI_Comm_rank(MPI_COMM_WORLD,&rank); if ( rank!=0 ) return; bool empty = map.created.empty() && map.destroyed.empty() && map.split.empty() && map.merge.empty() && map.merge_split.empty(); for (size_t i=0; i(map.created[i])); for (size_t i=0; i(map.destroyed[i])); for (size_t i=0; i(map.src_dst[i].first), static_cast(map.src_dst[i].second)); } for (size_t i=0; i(map.split[i].first), static_cast(map.split[i].second[0])); for (size_t j=1; j(map.split[i].second[j])); } for (size_t i=0; i(map.merge[i].first[0])); for (size_t j=1; j(map.merge[i].first[j])); fprintf(fid,"-%lli",static_cast(map.merge[i].second)); } for (size_t i=0; i(map.merge_split[i].first[0])); for (size_t j=1; j(map.merge_split[i].first[j])); fprintf(fid,"-%lli",static_cast(map.merge_split[i].second[0])); for (size_t j=1; j(map.merge_split[i].second[j])); } fprintf(fid,"\n"); fclose(fid); }