libMesh::Parallel::Sort< KeyType, IdxType > Class Template Reference

Object for performing parallel sorts using MPI. More...

#include <parallel_sort.h>

Inheritance diagram for libMesh::Parallel::Sort< KeyType, IdxType >:

Public Member Functions

 Sort (const Parallel::Communicator &comm, std::vector< KeyType > &d)
 
void sort ()
 
const std::vector< KeyType > & bin ()
 
const Parallel::Communicatorcomm () const
 
processor_id_type n_processors () const
 
processor_id_type processor_id () const
 

Protected Attributes

const Parallel::Communicator_communicator
 

Private Member Functions

void binsort ()
 
void communicate_bins ()
 
void sort_local_bin ()
 
template<>
void binsort ()
 

Private Attributes

const processor_id_type _n_procs
 
const processor_id_type _proc_id
 
bool _bin_is_sorted
 
std::vector< KeyType > & _data
 
std::vector< IdxType > _local_bin_sizes
 
std::vector< KeyType > _my_bin
 

Detailed Description

template<typename KeyType, typename IdxType = unsigned int>
class libMesh::Parallel::Sort< KeyType, IdxType >

Object for performing parallel sorts using MPI.

The parallel sorting method is templated on the type of data which is to be sorted. It may later be templated on other things if we are ambitious. This class knows about MPI, and knows how many processors there are. It is responsible for transmitting data between the processors and ensuring that the data is properly sorted between all the processors. We assume that a Sort is instantiated on all processors.

Author
Benjamin S. Kirk
John W. Peterson
Date
2007

Definition at line 53 of file parallel_sort.h.

Constructor & Destructor Documentation

◆ Sort()

template<typename KeyType , typename IdxType >
libMesh::Parallel::Sort< KeyType, IdxType >::Sort ( const Parallel::Communicator comm,
std::vector< KeyType > &  d 
)

Constructor takes the number of processors, the processor id, and a reference to a vector of data to be sorted. This vector is sorted by the constructor, therefore, construction of a Sort object takes O(n log n) time, where n is the length of the vector.

Definition at line 43 of file parallel_sort.C.

References libMesh::Parallel::Sort< KeyType, IdxType >::_data, libMesh::Parallel::Sort< KeyType, IdxType >::_local_bin_sizes, and libMesh::Parallel::Sort< KeyType, IdxType >::_n_procs.

44  :
45  ParallelObject(comm_in),
46  _n_procs(cast_int<processor_id_type>(comm_in.size())),
47  _proc_id(cast_int<processor_id_type>(comm_in.rank())),
48  _bin_is_sorted(false),
49  _data(d)
50 {
51  std::sort(_data.begin(), _data.end());
52 
53  // Allocate storage
54  _local_bin_sizes.resize(_n_procs);
55 }
ParallelObject(const Parallel::Communicator &comm_in)
const processor_id_type _proc_id
Definition: parallel_sort.h:93
std::vector< KeyType > & _data
std::vector< IdxType > _local_bin_sizes
const processor_id_type _n_procs
Definition: parallel_sort.h:88

Member Function Documentation

◆ bin()

template<typename KeyType , typename IdxType >
const std::vector< KeyType > & libMesh::Parallel::Sort< KeyType, IdxType >::bin ( )

Return a constant reference to _my_bin. This allows us to do things like check if sorting was successful by printing _my_bin.

Definition at line 260 of file parallel_sort.C.

References libMesh::out.

Referenced by libMesh::MeshCommunication::assign_global_indices(), and libMesh::MeshCommunication::find_global_indices().

261 {
262  if (!_bin_is_sorted)
263  {
264  libMesh::out << "Warning! Bin is not yet sorted!" << std::endl;
265  }
266 
267  return _my_bin;
268 }
std::vector< KeyType > _my_bin
OStreamProxy out(std::cout)

◆ binsort() [1/2]

template<typename KeyType , typename IdxType >
void libMesh::Parallel::Sort< KeyType, IdxType >::binsort ( )
private

Sorts the local data into bins across all processors. Right now it constructs a BenSorter<KeyType> object. In the future this could be a template parameter.

Definition at line 98 of file parallel_sort.C.

References libMesh::Parallel::BinSorter< KeyType, IdxType >::binsort(), and libMesh::Parallel::BinSorter< KeyType, IdxType >::sizeof_bin().

99 {
100  // Find the global min and max from all the
101  // processors.
102  std::vector<KeyType> global_min_max(2);
103 
104  // Insert the local min and max for this processor
105  global_min_max[0] = -_data.front();
106  global_min_max[1] = _data.back();
107 
108  // Communicate to determine the global
109  // min and max for all processors.
110  this->comm().max(global_min_max);
111 
112  // Multiply the min by -1 to obtain the true min
113  global_min_max[0] *= -1;
114 
115  // Bin-Sort based on the global min and max
116  Parallel::BinSorter<KeyType> bs(this->comm(), _data);
117  bs.binsort(_n_procs, global_min_max[1], global_min_max[0]);
118 
119  // Now save the local bin sizes in a vector so
120  // we don't have to keep around the BinSorter.
121  for (processor_id_type i=0; i<_n_procs; ++i)
122  _local_bin_sizes[i] = bs.sizeof_bin(i);
123 }
uint8_t processor_id_type
Definition: id_types.h:99
const Parallel::Communicator & comm() const
std::vector< KeyType > & _data
std::vector< IdxType > _local_bin_sizes
const processor_id_type _n_procs
Definition: parallel_sort.h:88

◆ binsort() [2/2]

template<>
void libMesh::Parallel::Sort< Parallel::DofObjectKey, unsigned int >::binsort ( )
private

Definition at line 132 of file parallel_sort.C.

References libMesh::Parallel::BinSorter< KeyType, IdxType >::binsort(), dofobjectkey_max_op(), dofobjectkey_min_op(), std::max(), and libMesh::Parallel::BinSorter< KeyType, IdxType >::sizeof_bin().

133 {
134  // Find the global min and max from all the
135  // processors. Do this using MPI_Allreduce.
137  local_min, local_max,
138  global_min, global_max;
139 
140  if (_data.empty())
141  {
142 #ifdef LIBMESH_ENABLE_UNIQUE_ID
143  local_min.first.rack0 = local_min.first.rack1 = local_min.first.rack2 = static_cast<Hilbert::inttype>(-1);
144  local_min.second = std::numeric_limits<unique_id_type>::max();
145  local_max.first.rack0 = local_max.first.rack1 = local_max.first.rack2 = 0;
146  local_max.second = 0;
147 #else
148  local_min.rack0 = local_min.rack1 = local_min.rack2 = static_cast<Hilbert::inttype>(-1);
149  local_max.rack0 = local_max.rack1 = local_max.rack2 = 0;
150 #endif
151  }
152  else
153  {
154  local_min = _data.front();
155  local_max = _data.back();
156  }
157 
158  MPI_Op hilbert_max, hilbert_min;
159 
160  MPI_Op_create ((MPI_User_function*)dofobjectkey_max_op, true, &hilbert_max);
161  MPI_Op_create ((MPI_User_function*)dofobjectkey_min_op, true, &hilbert_min);
162 
163  // Communicate to determine the global
164  // min and max for all processors.
165  MPI_Allreduce(&local_min,
166  &global_min,
167  1,
168  Parallel::StandardType<Parallel::DofObjectKey>(&local_min),
169  hilbert_min,
170  this->comm().get());
171 
172  MPI_Allreduce(&local_max,
173  &global_max,
174  1,
175  Parallel::StandardType<Parallel::DofObjectKey>(&local_max),
176  hilbert_max,
177  this->comm().get());
178 
179  MPI_Op_free (&hilbert_max);
180  MPI_Op_free (&hilbert_min);
181 
182  // Bin-Sort based on the global min and max
183  Parallel::BinSorter<Parallel::DofObjectKey> bs(this->comm(),_data);
184  bs.binsort(_n_procs, global_max, global_min);
185 
186  // Now save the local bin sizes in a vector so
187  // we don't have to keep around the BinSorter.
188  for (processor_id_type i=0; i<_n_procs; ++i)
189  _local_bin_sizes[i] = bs.sizeof_bin(i);
190 }
void dofobjectkey_min_op(libMesh::Parallel::DofObjectKey *in, libMesh::Parallel::DofObjectKey *inout, int *len, void *)
std::pair< Hilbert::HilbertIndices, unique_id_type > DofObjectKey
uint8_t processor_id_type
Definition: id_types.h:99
const Parallel::Communicator & comm() const
std::vector< KeyType > & _data
std::vector< IdxType > _local_bin_sizes
long double max(long double a, double b)
void dofobjectkey_max_op(libMesh::Parallel::DofObjectKey *in, libMesh::Parallel::DofObjectKey *inout, int *len, void *)
const processor_id_type _n_procs
Definition: parallel_sort.h:88

◆ comm()

const Parallel::Communicator& libMesh::ParallelObject::comm ( ) const
inlineinherited
Returns
A reference to the Parallel::Communicator object used by this mesh.

Definition at line 89 of file parallel_object.h.

References libMesh::ParallelObject::_communicator.

Referenced by libMesh::__libmesh_petsc_diff_solver_jacobian(), libMesh::__libmesh_petsc_diff_solver_monitor(), libMesh::__libmesh_petsc_diff_solver_residual(), libMesh::__libmesh_tao_equality_constraints(), libMesh::__libmesh_tao_equality_constraints_jacobian(), libMesh::__libmesh_tao_gradient(), libMesh::__libmesh_tao_hessian(), libMesh::__libmesh_tao_inequality_constraints(), libMesh::__libmesh_tao_inequality_constraints_jacobian(), libMesh::__libmesh_tao_objective(), libMesh::MeshRefinement::_coarsen_elements(), libMesh::ExactSolution::_compute_error(), libMesh::UniformRefinementEstimator::_estimate_error(), libMesh::BoundaryInfo::_find_id_maps(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_get_diagonal(), libMesh::SlepcEigenSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult(), libMesh::PetscLinearSolver< T >::_petsc_shell_matrix_mult_add(), libMesh::EquationSystems::_read_impl(), libMesh::MeshRefinement::_refine_elements(), libMesh::MeshRefinement::_smooth_flags(), libMesh::PetscDMWrapper::add_dofs_helper(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::ImplicitSystem::add_matrix(), libMesh::System::add_vector(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assemble_qoi(), libMesh::MeshCommunication::assign_global_indices(), libMesh::DofMap::attach_matrix(), libMesh::MeshTools::Generation::build_extrusion(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::EquationSystems::build_parallel_solution_vector(), libMesh::PetscDMWrapper::build_section(), libMesh::PetscDMWrapper::build_sf(), libMesh::MeshBase::cache_elem_dims(), libMesh::System::calculate_norm(), libMesh::DofMap::check_dirichlet_bcid_consistency(), libMesh::PetscDMWrapper::check_section_n_dofs(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Problem_Interface::computeF(), libMesh::Problem_Interface::computeJacobian(), libMesh::Problem_Interface::computePreconditioner(), libMesh::ExodusII_IO::copy_elemental_solution(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::MeshTools::create_bounding_box(), libMesh::MeshTools::create_nodal_bounding_box(), libMesh::MeshRefinement::create_parent_error_vector(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::MeshTools::create_subdomain_bounding_box(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), DMlibMeshFunction(), DMlibMeshJacobian(), DMlibMeshSetSystem_libMesh(), DMVariableBounds_libMesh(), libMesh::MeshRefinement::eliminate_unrefined_patches(), libMesh::EpetraVector< T >::EpetraVector(), libMesh::WeightedPatchRecoveryErrorEstimator::estimate_error(), libMesh::PatchRecoveryErrorEstimator::estimate_error(), libMesh::JumpErrorEstimator::estimate_error(), libMesh::AdjointRefinementEstimator::estimate_error(), libMesh::ExactErrorEstimator::estimate_error(), libMesh::MeshRefinement::flag_elements_by_elem_fraction(), libMesh::MeshRefinement::flag_elements_by_error_fraction(), libMesh::MeshRefinement::flag_elements_by_nelem_target(), libMesh::CondensedEigenSystem::get_eigenpair(), libMesh::DofMap::get_info(), libMesh::ImplicitSystem::get_linear_solver(), libMesh::LocationMap< T >::init(), libMesh::TimeSolver::init(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::EigenSystem::init_data(), libMesh::EigenSystem::init_matrices(), libMesh::OptimizationSystem::initialize_equality_constraints_storage(), libMesh::OptimizationSystem::initialize_inequality_constraints_storage(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_flags(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_p_levels(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::MeshTools::libmesh_assert_valid_unique_ids(), libMesh::libmesh_petsc_snes_fd_residual(), libMesh::libmesh_petsc_snes_jacobian(), libMesh::libmesh_petsc_snes_mffd_residual(), libMesh::libmesh_petsc_snes_postcheck(), libMesh::libmesh_petsc_snes_residual(), libMesh::libmesh_petsc_snes_residual_helper(), libMesh::MeshRefinement::limit_level_mismatch_at_edge(), libMesh::MeshRefinement::limit_level_mismatch_at_node(), libMesh::MeshRefinement::limit_overrefined_boundary(), libMesh::MeshRefinement::limit_underrefined_boundary(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshCommunication::make_elems_parallel_consistent(), libMesh::MeshRefinement::make_flags_parallel_consistent(), libMesh::MeshCommunication::make_new_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_new_nodes_parallel_consistent(), libMesh::MeshCommunication::make_node_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_proc_ids_parallel_consistent(), libMesh::MeshCommunication::make_node_unique_ids_parallel_consistent(), libMesh::MeshCommunication::make_nodes_parallel_consistent(), libMesh::MeshCommunication::make_p_levels_parallel_consistent(), libMesh::MeshRefinement::make_refinement_compatible(), libMesh::FEMSystem::mesh_position_set(), libMesh::DistributedMesh::n_active_elem(), libMesh::MeshTools::n_active_levels(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::CondensedEigenSystem::n_global_non_condensed_dofs(), libMesh::MeshTools::n_levels(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::MeshTools::n_p_levels(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::DistributedMesh::parallel_max_elem_id(), libMesh::DistributedMesh::parallel_max_node_id(), libMesh::ReplicatedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_max_unique_id(), libMesh::DistributedMesh::parallel_n_elem(), libMesh::DistributedMesh::parallel_n_nodes(), libMesh::SparsityPattern::Build::parallel_sync(), libMesh::MeshTools::paranoid_n_levels(), libMesh::petsc_auto_fieldsplit(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::MeshBase::prepare_for_use(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::MeshBase::recalculate_n_partitions(), libMesh::MeshRefinement::refine_and_coarsen_elements(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::PetscDiffSolver::setup_petsc_data(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::split_mesh(), libMesh::MeshBase::subdomain_ids(), libMesh::BoundaryInfo::sync(), libMesh::MeshRefinement::test_level_one(), libMesh::MeshRefinement::test_unflagged(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

90  { return _communicator; }
const Parallel::Communicator & _communicator

◆ communicate_bins()

template<typename KeyType , typename IdxType >
void libMesh::Parallel::Sort< KeyType, IdxType >::communicate_bins ( )
private

Communicates the bins from each processor to the appropriate processor. By the time this function is finished, each processor will hold only its own bin(s).

Definition at line 196 of file parallel_sort.C.

References end, and libMesh::Parallel::push_parallel_vector_data().

197 {
198 #ifdef LIBMESH_HAVE_MPI
199  // Find each section of our data to send
200  IdxType local_offset = 0;
201  std::map<processor_id_type, std::vector<KeyType> > pushed_keys, received_keys;
202 
203  for (processor_id_type i=0; i != _n_procs; ++i)
204  {
205  IdxType next_offset = local_offset + _local_bin_sizes[i];
206  if (_local_bin_sizes[i])
207  {
208  auto begin = _data.begin() + local_offset;
209  auto end = _data.begin() + next_offset;
210  pushed_keys[i].assign(begin, end);
211  }
212 
213  local_offset = next_offset;
214  }
215 
216  auto keys_action_functor =
217  [& received_keys]
218  (processor_id_type pid,
219  const std::vector<KeyType> & keys)
220  {
221  received_keys[pid] = keys;
222  };
223 
225  (this->comm(), pushed_keys, keys_action_functor);
226 
227  std::size_t my_bin_size = 0;
228  for (auto & p : received_keys)
229  my_bin_size += p.second.size();
230 
231  _my_bin.clear();
232  _my_bin.reserve(my_bin_size);
233 
234  for (auto & p : received_keys)
235  _my_bin.insert(_my_bin.end(), p.second.begin(), p.second.end());
236 
237 #ifdef DEBUG
238  std::vector<IdxType> global_bin_sizes = _local_bin_sizes;
239 
240  this->comm().sum(global_bin_sizes);
241 
242  libmesh_assert_equal_to
243  (global_bin_sizes[this->processor_id()], _my_bin.size());
244 #endif
245 
246 #endif // LIBMESH_HAVE_MPI
247 }
uint8_t processor_id_type
Definition: id_types.h:99
const Parallel::Communicator & comm() const
std::vector< KeyType > & _data
IterBase * end
std::vector< IdxType > _local_bin_sizes
const processor_id_type _n_procs
Definition: parallel_sort.h:88
std::vector< KeyType > _my_bin
void push_parallel_vector_data(const Communicator &comm, const MapToVectors &data, RequestContainer &reqs, ActionFunctor &act_on_data)
processor_id_type processor_id() const

◆ n_processors()

processor_id_type libMesh::ParallelObject::n_processors ( ) const
inlineinherited
Returns
The number of processors in the group.

Definition at line 95 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::size().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::DistributedMesh::add_node(), libMesh::LaplaceMeshSmoother::allgather_graph(), libMesh::FEMSystem::assembly(), libMesh::AztecLinearSolver< T >::AztecLinearSolver(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::UnstructuredMesh::create_pid_mesh(), libMesh::MeshTools::create_processor_bounding_box(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshBase::get_info(), libMesh::SystemSubsetBySubdomain::init(), libMesh::PetscDMWrapper::init_and_attach_petscdm(), libMesh::Nemesis_IO_Helper::initialize(), libMesh::DistributedMesh::insert_elem(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_new_node_procids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_topology_consistent_procids< Node >(), libMesh::MeshTools::libmesh_assert_valid_boundary_ids(), libMesh::MeshTools::libmesh_assert_valid_dof_ids(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::MeshTools::libmesh_assert_valid_refinement_flags(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_elem_on_proc(), libMesh::MeshBase::n_elem_on_proc(), libMesh::MeshBase::n_nodes_on_proc(), libMesh::MeshBase::partition(), libMesh::PetscLinearSolver< T >::PetscLinearSolver(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::CheckpointIO::read(), libMesh::CheckpointIO::read_connectivity(), libMesh::XdrIO::read_header(), libMesh::CheckpointIO::read_nodes(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::System::read_serialized_vector(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::GMVIO::write_binary(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::XdrIO::write_serialized_nodes(), and libMesh::XdrIO::write_serialized_nodesets().

96  { return cast_int<processor_id_type>(_communicator.size()); }
processor_id_type size() const
Definition: communicator.h:175
const Parallel::Communicator & _communicator

◆ processor_id()

processor_id_type libMesh::ParallelObject::processor_id ( ) const
inlineinherited
Returns
The rank of this processor in the group.

Definition at line 101 of file parallel_object.h.

References libMesh::ParallelObject::_communicator, and libMesh::Parallel::Communicator::rank().

Referenced by libMesh::BoundaryInfo::_find_id_maps(), libMesh::EquationSystems::_read_impl(), libMesh::PetscDMWrapper::add_dofs_to_section(), libMesh::DistributedMesh::add_elem(), libMesh::BoundaryInfo::add_elements(), libMesh::DofMap::add_neighbors_to_send_list(), libMesh::DistributedMesh::add_node(), libMesh::UnstructuredMesh::all_second_order(), libMesh::MeshTools::Modification::all_tri(), libMesh::FEMSystem::assembly(), libMesh::EquationSystems::build_discontinuous_solution_vector(), libMesh::Nemesis_IO_Helper::build_element_and_node_maps(), libMesh::InfElemBuilder::build_inf_elem(), libMesh::BoundaryInfo::build_node_list_from_side_list(), libMesh::EquationSystems::build_parallel_elemental_solution_vector(), libMesh::DistributedMesh::clear(), libMesh::ExodusII_IO_Helper::close(), libMesh::Nemesis_IO_Helper::compute_border_node_ids(), libMesh::Nemesis_IO_Helper::compute_communication_map_parameters(), libMesh::Nemesis_IO_Helper::compute_internal_and_border_elems_and_internal_nodes(), libMesh::Nemesis_IO_Helper::compute_node_communication_maps(), libMesh::Nemesis_IO_Helper::compute_num_global_elem_blocks(), libMesh::Nemesis_IO_Helper::compute_num_global_nodesets(), libMesh::Nemesis_IO_Helper::compute_num_global_sidesets(), libMesh::Nemesis_IO_Helper::construct_nemesis_filename(), libMesh::MeshTools::correct_node_proc_ids(), libMesh::ExodusII_IO_Helper::create(), libMesh::DistributedMesh::delete_elem(), libMesh::DistributedMesh::delete_node(), libMesh::MeshCommunication::delete_remote_elements(), libMesh::DofMap::distribute_dofs(), libMesh::DofMap::distribute_local_dofs_node_major(), libMesh::DofMap::distribute_local_dofs_var_major(), libMesh::DistributedMesh::DistributedMesh(), libMesh::DofMap::end_dof(), libMesh::DofMap::end_old_dof(), libMesh::EnsightIO::EnsightIO(), libMesh::MeshFunction::find_element(), libMesh::MeshFunction::find_elements(), libMesh::UnstructuredMesh::find_neighbors(), libMesh::DofMap::first_dof(), libMesh::DofMap::first_old_dof(), libMesh::Nemesis_IO_Helper::get_cmap_params(), libMesh::Nemesis_IO_Helper::get_eb_info_global(), libMesh::Nemesis_IO_Helper::get_elem_cmap(), libMesh::Nemesis_IO_Helper::get_elem_map(), libMesh::MeshBase::get_info(), libMesh::DofMap::get_info(), libMesh::Nemesis_IO_Helper::get_init_global(), libMesh::Nemesis_IO_Helper::get_init_info(), libMesh::Nemesis_IO_Helper::get_loadbal_param(), libMesh::Nemesis_IO_Helper::get_node_cmap(), libMesh::Nemesis_IO_Helper::get_node_map(), libMesh::Nemesis_IO_Helper::get_ns_param_global(), libMesh::Nemesis_IO_Helper::get_ss_param_global(), libMesh::SparsityPattern::Build::handle_vi_vj(), libMesh::SystemSubsetBySubdomain::init(), libMesh::ExodusII_IO_Helper::initialize(), libMesh::ExodusII_IO_Helper::initialize_element_variables(), libMesh::ExodusII_IO_Helper::initialize_global_variables(), libMesh::ExodusII_IO_Helper::initialize_nodal_variables(), libMesh::DistributedMesh::insert_elem(), libMesh::DofMap::is_evaluable(), libMesh::SparsityPattern::Build::join(), libMesh::DofMap::last_dof(), libMesh::MeshTools::libmesh_assert_consistent_distributed(), libMesh::MeshTools::libmesh_assert_consistent_distributed_nodes(), libMesh::MeshTools::libmesh_assert_contiguous_dof_ids(), libMesh::MeshTools::libmesh_assert_parallel_consistent_procids< Elem >(), libMesh::MeshTools::libmesh_assert_valid_neighbors(), libMesh::DistributedMesh::libmesh_assert_valid_parallel_object_ids(), libMesh::DofMap::local_variable_indices(), libMesh::MeshRefinement::make_coarsening_compatible(), libMesh::MeshBase::n_active_local_elem(), libMesh::BoundaryInfo::n_boundary_conds(), libMesh::BoundaryInfo::n_edge_conds(), libMesh::DofMap::n_local_dofs(), libMesh::System::n_local_dofs(), libMesh::MeshBase::n_local_elem(), libMesh::MeshBase::n_local_nodes(), libMesh::BoundaryInfo::n_nodeset_conds(), libMesh::BoundaryInfo::n_shellface_conds(), libMesh::SparsityPattern::Build::operator()(), libMesh::DistributedMesh::own_node(), libMesh::System::point_gradient(), libMesh::System::point_hessian(), libMesh::System::point_value(), libMesh::Nemesis_IO_Helper::put_cmap_params(), libMesh::Nemesis_IO_Helper::put_elem_cmap(), libMesh::Nemesis_IO_Helper::put_elem_map(), libMesh::Nemesis_IO_Helper::put_loadbal_param(), libMesh::Nemesis_IO_Helper::put_node_cmap(), libMesh::Nemesis_IO_Helper::put_node_map(), libMesh::NameBasedIO::read(), libMesh::Nemesis_IO::read(), libMesh::XdrIO::read(), libMesh::CheckpointIO::read(), libMesh::ExodusII_IO_Helper::read_elem_num_map(), libMesh::ExodusII_IO_Helper::read_global_values(), libMesh::CheckpointIO::read_header(), libMesh::XdrIO::read_header(), libMesh::System::read_header(), libMesh::System::read_legacy_data(), libMesh::ExodusII_IO_Helper::read_node_num_map(), libMesh::System::read_parallel_data(), libMesh::System::read_SCALAR_dofs(), libMesh::XdrIO::read_serialized_bc_names(), libMesh::XdrIO::read_serialized_bcs_helper(), libMesh::System::read_serialized_blocked_dof_objects(), libMesh::XdrIO::read_serialized_connectivity(), libMesh::System::read_serialized_data(), libMesh::XdrIO::read_serialized_nodes(), libMesh::XdrIO::read_serialized_nodesets(), libMesh::XdrIO::read_serialized_subdomain_names(), libMesh::System::read_serialized_vector(), libMesh::System::read_serialized_vectors(), libMesh::DistributedMesh::renumber_dof_objects(), libMesh::CheckpointIO::select_split_config(), libMesh::DofMap::set_nonlocal_dof_objects(), libMesh::PetscDMWrapper::set_point_range_in_section(), libMesh::LaplaceMeshSmoother::smooth(), libMesh::MeshTools::total_weight(), libMesh::MeshRefinement::uniformly_coarsen(), libMesh::Parallel::Packing< T >::unpack(), libMesh::DistributedMesh::update_parallel_id_counts(), libMesh::NameBasedIO::write(), libMesh::XdrIO::write(), libMesh::CheckpointIO::write(), libMesh::EquationSystems::write(), libMesh::GMVIO::write_discontinuous_gmv(), libMesh::ExodusII_IO::write_element_data(), libMesh::ExodusII_IO_Helper::write_element_values(), libMesh::ExodusII_IO_Helper::write_elements(), libMesh::ExodusII_IO::write_global_data(), libMesh::ExodusII_IO_Helper::write_global_values(), libMesh::System::write_header(), libMesh::ExodusII_IO::write_information_records(), libMesh::ExodusII_IO_Helper::write_information_records(), libMesh::ExodusII_IO_Helper::write_nodal_coordinates(), libMesh::UCDIO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data(), libMesh::ExodusII_IO::write_nodal_data_discontinuous(), libMesh::ExodusII_IO_Helper::write_nodal_values(), libMesh::Nemesis_IO_Helper::write_nodesets(), libMesh::ExodusII_IO_Helper::write_nodesets(), libMesh::System::write_parallel_data(), libMesh::System::write_SCALAR_dofs(), libMesh::XdrIO::write_serialized_bc_names(), libMesh::XdrIO::write_serialized_bcs_helper(), libMesh::System::write_serialized_blocked_dof_objects(), libMesh::XdrIO::write_serialized_connectivity(), libMesh::System::write_serialized_data(), libMesh::XdrIO::write_serialized_nodes(), libMesh::XdrIO::write_serialized_nodesets(), libMesh::XdrIO::write_serialized_subdomain_names(), libMesh::System::write_serialized_vector(), libMesh::System::write_serialized_vectors(), libMesh::Nemesis_IO_Helper::write_sidesets(), libMesh::ExodusII_IO_Helper::write_sidesets(), libMesh::ExodusII_IO::write_timestep(), libMesh::ExodusII_IO_Helper::write_timestep(), and libMesh::ExodusII_IO::write_timestep_discontinuous().

102  { return cast_int<processor_id_type>(_communicator.rank()); }
const Parallel::Communicator & _communicator
processor_id_type rank() const
Definition: communicator.h:173

◆ sort()

template<typename KeyType , typename IdxType >
void libMesh::Parallel::Sort< KeyType, IdxType >::sort ( )

This is the only method which needs to be called by the user. Its only responsibility is to call three private methods in the correct order.

Definition at line 60 of file parallel_sort.C.

Referenced by libMesh::MeshCommunication::assign_global_indices(), and libMesh::MeshCommunication::find_global_indices().

61 {
62  // Find the global data size. The sorting
63  // algorithms assume they have a range to
64  // work with, so catch the degenerate cases here
65  IdxType global_data_size = cast_int<IdxType>(_data.size());
66 
67  this->comm().sum (global_data_size);
68 
69  if (global_data_size < 2)
70  {
71  // the entire global range is either empty
72  // or contains only one element
73  _my_bin = _data;
74 
75  this->comm().allgather (static_cast<IdxType>(_my_bin.size()),
77  }
78  else
79  {
80  if (this->n_processors() > 1)
81  {
82  this->binsort();
83  this->communicate_bins();
84  }
85  else
86  _my_bin = _data;
87 
88  this->sort_local_bin();
89  }
90 
91  // Set sorted flag to true
92  _bin_is_sorted = true;
93 }
void allgather(const T &send, std::vector< T, A > &recv) const
const Parallel::Communicator & comm() const
std::vector< KeyType > & _data
std::vector< IdxType > _local_bin_sizes
processor_id_type n_processors() const
std::vector< KeyType > _my_bin

◆ sort_local_bin()

template<typename KeyType , typename IdxType >
void libMesh::Parallel::Sort< KeyType, IdxType >::sort_local_bin ( )
private

After all the bins have been communicated, we can sort our local bin. This is nothing more than a call to std::sort

Definition at line 252 of file parallel_sort.C.

253 {
254  std::sort(_my_bin.begin(), _my_bin.end());
255 }
std::vector< KeyType > _my_bin

Member Data Documentation

◆ _bin_is_sorted

template<typename KeyType, typename IdxType = unsigned int>
bool libMesh::Parallel::Sort< KeyType, IdxType >::_bin_is_sorted
private

Flag which lets you know if sorting is complete

Definition at line 98 of file parallel_sort.h.

◆ _communicator

◆ _data

template<typename KeyType, typename IdxType = unsigned int>
std::vector<KeyType>& libMesh::Parallel::Sort< KeyType, IdxType >::_data
private

The raw, unsorted data which will need to be sorted (in parallel) across all processors.

Definition at line 105 of file parallel_sort.h.

Referenced by libMesh::Parallel::Sort< KeyType, IdxType >::Sort().

◆ _local_bin_sizes

template<typename KeyType, typename IdxType = unsigned int>
std::vector<IdxType> libMesh::Parallel::Sort< KeyType, IdxType >::_local_bin_sizes
private

Vector which holds the size of each bin on this processor. It has size equal to _n_procs.

Definition at line 112 of file parallel_sort.h.

Referenced by libMesh::Parallel::Sort< KeyType, IdxType >::Sort().

◆ _my_bin

template<typename KeyType, typename IdxType = unsigned int>
std::vector<KeyType> libMesh::Parallel::Sort< KeyType, IdxType >::_my_bin
private

The bin which will eventually be held by this processor. It may be shorter or longer than _data. It will be dynamically resized when it is needed.

Definition at line 120 of file parallel_sort.h.

◆ _n_procs

template<typename KeyType, typename IdxType = unsigned int>
const processor_id_type libMesh::Parallel::Sort< KeyType, IdxType >::_n_procs
private

The number of processors to work with.

Definition at line 88 of file parallel_sort.h.

Referenced by libMesh::Parallel::Sort< KeyType, IdxType >::Sort().

◆ _proc_id

template<typename KeyType, typename IdxType = unsigned int>
const processor_id_type libMesh::Parallel::Sort< KeyType, IdxType >::_proc_id
private

The identity of this processor.

Definition at line 93 of file parallel_sort.h.


The documentation for this class was generated from the following files: