Loading [MathJax]/jax/output/HTML-CSS/config.js
OGS
BaseLib::MPI Namespace Reference

Classes

struct  Mpi
 
struct  Setup
 

Functions

template<typename T >
constexpr MPI_Datatype mpiType ()
 
template<typename T >
static std::vector< T > allgather (T const &value, Mpi const &mpi)
 
template<typename T >
static std::vector< T > allgather (std::vector< T > const &vector, Mpi const &mpi)
 
template<typename T >
static T allreduce (T const &value, MPI_Op const &mpi_op, Mpi const &mpi)
 
template<typename T >
static std::vector< T > allreduce (std::vector< T > const &vector, MPI_Op const &mpi_op, Mpi const &mpi)
 
template<typename T >
static void allreduceInplace (std::vector< T > &vector, MPI_Op const &mpi_op, Mpi const &mpi)
 
template<typename T >
static std::vector< int > allgatherv (std::span< T > const send_buffer, std::vector< std::remove_const_t< T > > &receive_buffer, Mpi const &mpi)
 
static bool anyOf (bool const val, Mpi const &mpi=Mpi{MPI_COMM_WORLD})
 

Function Documentation

◆ allgather() [1/2]

template<typename T >
static std::vector< T > BaseLib::MPI::allgather ( std::vector< T > const & vector,
Mpi const & mpi )
static

Definition at line 111 of file MPI.h.

112{
113 std::size_t const size = vector.size();
114 // Flat in memory over all ranks;
115 std::vector<T> result(mpi.size * size);
116
117 MPI_Allgather(vector.data(), size, mpiType<T>(), result.data(), size,
118 mpiType<T>(), mpi.communicator);
119
120 return result;
121}

References BaseLib::MPI::Mpi::communicator, mpiType(), and BaseLib::MPI::Mpi::size.

◆ allgather() [2/2]

template<typename T >
static std::vector< T > BaseLib::MPI::allgather ( T const & value,
Mpi const & mpi )
static

Definition at line 100 of file MPI.h.

101{
102 std::vector<T> result(mpi.size);
103
104 MPI_Allgather(&value, 1, mpiType<T>(), result.data(), 1, mpiType<T>(),
105 mpi.communicator);
106
107 return result;
108}

References BaseLib::MPI::Mpi::communicator, mpiType(), and BaseLib::MPI::Mpi::size.

Referenced by allgatherv(), MeshLib::computeNumberOfRegularBaseNodesAtRank(), MeshLib::computeNumberOfRegularHigherOrderNodesAtRank(), MeshLib::computeRegularBaseNodeGlobalNodeIDsOfSubDomainPartition(), MeshLib::IO::getPartitionInfo(), and MeshLib::IO::NodePartitionedMeshReader::newMesh().

◆ allgatherv()

template<typename T >
static std::vector< int > BaseLib::MPI::allgatherv ( std::span< T > const send_buffer,
std::vector< std::remove_const_t< T > > & receive_buffer,
Mpi const & mpi )
static

Allgather for variable data. Returns offsets in the receive buffer. The receive buffer is resized to accommodate the gathered data.

Definition at line 160 of file MPI.h.

164{
165 // Determine the number of elements to send
166 int const size = static_cast<int>(send_buffer.size());
167
168 // Gather sizes from all ranks
169 std::vector<int> const sizes = allgather(size, mpi);
170
171 // Compute offsets based on counts
172 std::vector<int> const offsets = BaseLib::sizesToOffsets(sizes);
173
174 // Resize receive buffer to hold all gathered data
175 receive_buffer.resize(offsets.back());
176
177 MPI_Allgatherv(send_buffer.data(), size, mpiType<T>(),
178 receive_buffer.data(), sizes.data(), offsets.data(),
179 mpiType<T>(), mpi.communicator);
180
181 return offsets;
182}
static std::vector< T > allgather(T const &value, Mpi const &mpi)
Definition MPI.h:100
std::vector< ranges::range_value_t< R > > sizesToOffsets(R const &sizes)
Definition Algorithm.h:283

References allgather(), BaseLib::MPI::Mpi::communicator, mpiType(), and BaseLib::sizesToOffsets().

Referenced by MeshLib::computeGhostBaseNodeGlobalNodeIDsOfSubDomainPartition(), and MathLib::PETScVector::getGlobalVector().

◆ allreduce() [1/2]

template<typename T >
static std::vector< T > BaseLib::MPI::allreduce ( std::vector< T > const & vector,
MPI_Op const & mpi_op,
Mpi const & mpi )
static

Definition at line 133 of file MPI.h.

135{
136 std::size_t const size = vector.size();
137 std::vector<T> result(vector.size());
138
139 MPI_Allreduce(vector.data(), result.data(), size, mpiType<T>(), mpi_op,
140 mpi.communicator);
141 return result;
142}

References BaseLib::MPI::Mpi::communicator, and mpiType().

◆ allreduce() [2/2]

◆ allreduceInplace()

template<typename T >
static void BaseLib::MPI::allreduceInplace ( std::vector< T > & vector,
MPI_Op const & mpi_op,
Mpi const & mpi )
static

Definition at line 145 of file MPI.h.

148{
149 MPI_Allreduce(MPI_IN_PLACE,
150 vector.data(),
151 vector.size(),
152 mpiType<T>(),
153 mpi_op,
154 mpi.communicator);
155}

References BaseLib::MPI::Mpi::communicator, and mpiType().

Referenced by MeshLib::computeGhostBaseNodeGlobalNodeIDsOfSubDomainPartition().

◆ anyOf()

static bool BaseLib::MPI::anyOf ( bool const val,
Mpi const & mpi = Mpi{MPI_COMM_WORLD} )
inlinestatic

The reduction is implemented transparently for with and without MPI. In the latter case the input value is returned.

Definition at line 187 of file MPI.h.

190 {MPI_COMM_WORLD}
191#endif
192)
193{
194#ifdef USE_PETSC
195 return allreduce(val, MPI_LOR, mpi);
196#else
197 return val;
198#endif
199}

Referenced by ProcessLib::AssemblyMixin< Process >::assembleWithJacobian(), and NumLib::NonlinearSolver< NonlinearSolverTag::Newton >::solve().

◆ mpiType()

template<typename T >
MPI_Datatype BaseLib::MPI::mpiType ( )
constexpr

Definition at line 66 of file MPI.h.

67{
68 using U = std::remove_const_t<T>;
69 if constexpr (std::is_same_v<U, bool>)
70 {
71 return MPI_C_BOOL;
72 }
73 if constexpr (std::is_same_v<U, char>)
74 {
75 return MPI_CHAR;
76 }
77 if constexpr (std::is_same_v<U, double>)
78 {
79 return MPI_DOUBLE;
80 }
81 if constexpr (std::is_same_v<U, float>)
82 {
83 return MPI_FLOAT;
84 }
85 if constexpr (std::is_same_v<U, int>)
86 {
87 return MPI_INT;
88 }
89 if constexpr (std::is_same_v<U, std::size_t>)
90 {
91 return MPI_UNSIGNED_LONG;
92 }
93 if constexpr (std::is_same_v<U, unsigned int>)
94 {
95 return MPI_UNSIGNED;
96 }
97}

Referenced by allgather(), allgather(), allgatherv(), allreduce(), allreduce(), and allreduceInplace().