QMCPACK
Communicate Class Reference

Wrapping information on parallelism. More...

+ Inheritance diagram for Communicate:
+ Collaboration diagram for Communicate:

Public Member Functions

 Communicate ()
 constructor More...
 
 Communicate (const Communicate &in_comm, int nparts)
 constructor that splits in_comm More...
 
virtual ~Communicate ()
 destructor Call proper finalization of Communication library More...
 
 Communicate (const Communicate &)=delete
 disable constructor More...
 
void initialize (int argc, char **argv)
 
Communicate NodeComm () const
 provide a node/shared-memory communicator from current (parent) communicator More...
 
void finalize ()
 
void barrier () const
 
void abort () const
 
void barrier_and_abort (const std::string &msg) const
 
void set_world ()
 
mpi_comm_type getMPI () const
 return the Communicator ID (typically MPI_WORLD_COMM) More...
 
int rank () const
 return the rank More...
 
int size () const
 return the number of tasks More...
 
int getGroupID () const
 return the group id More...
 
int getNumGroups () const
 return the number of intra_comms which belong to the same group More...
 
void cleanupMessage (void *)
 
void setNodeID (int i)
 
void setNumNodes (int n)
 
void setName (const std::string &aname)
 
void setName (const char *aname, int alen)
 
const std::string & getName () const
 
bool isGroupLeader ()
 return true if the current MPI rank is the group lead More...
 
template<typename T >
void allreduce (T &)
 
template<typename T >
void reduce (T &)
 
template<typename T >
void reduce (T *restrict, T *restrict, int n)
 
template<typename T >
void reduce_in_place (T *restrict, int n)
 
template<typename T >
void bcast (T &)
 
template<typename T >
void bcast (T *restrict, int n)
 
template<typename T >
void send (int dest, int tag, T &)
 
template<typename T >
void gather (T &sb, T &rb, int dest=0)
 
template<typename T , typename IT >
void gatherv (T &sb, T &rb, IT &counts, IT &displ, int dest=0)
 
template<typename T >
void allgather (T &sb, T &rb, int count)
 
template<typename T , typename IT >
void allgatherv (T &sb, T &rb, IT &counts, IT &displ)
 
template<typename T >
void scatter (T &sb, T &rb, int dest=0)
 
template<typename T , typename IT >
void scatterv (T &sb, T &rb, IT &counts, IT &displ, int source=0)
 
template<typename T >
request irecv (int source, int tag, T &)
 
template<typename T >
request isend (int dest, int tag, T &)
 
template<typename T >
request irecv (int source, int tag, T *, int n)
 
template<typename T >
request isend (int dest, int tag, T *, int n)
 
template<typename T , typename IT >
void gatherv (T *sb, T *rb, int n, IT &counts, IT &displ, int dest=0)
 
template<typename T , typename TMPI , typename IT >
void gatherv_in_place (T *buf, TMPI &datatype, IT &counts, IT &displ, int dest=0)
 
template<typename T >
void allgather (T *sb, T *rb, int count)
 
template<typename T >
void gsum (T &)
 
CommunicategetGroupLeaderComm ()
 
template<>
void allreduce (int &g)
 
template<>
void allreduce (long &g)
 
template<>
void allreduce (unsigned long &g)
 
template<>
void allreduce (float &g)
 
template<>
void allreduce (double &g)
 
template<>
void allreduce (qmcplusplus::TinyVector< float, OHMMS_DIM > &g)
 
template<>
void allreduce (qmcplusplus::TinyVector< double, OHMMS_DIM > &g)
 
template<>
void allreduce (qmcplusplus::TinyVector< int, OHMMS_DIM > &g)
 
template<>
void allreduce (std::vector< int > &g)
 
template<>
void allreduce (std::vector< long > &g)
 
template<>
void allreduce (std::vector< unsigned long > &g)
 
template<>
void allreduce (std::vector< float > &g)
 
template<>
void allreduce (std::vector< double > &g)
 
template<>
void allreduce (std::vector< std::complex< float >> &g)
 
template<>
void allreduce (std::vector< std::complex< double >> &g)
 
template<>
void allreduce (PooledData< float > &g)
 
template<>
void allreduce (PooledData< double > &g)
 
template<>
void allreduce (qmcplusplus::Matrix< float > &g)
 
template<>
void allreduce (qmcplusplus::Matrix< double > &g)
 
template<>
void reduce (std::vector< float > &g)
 
template<>
void reduce (std::vector< double > &g)
 
template<>
void reduce (std::vector< int > &g)
 
template<>
void reduce (std::vector< long > &g)
 
template<>
void reduce (int *restrict g, int *restrict res, int n)
 
template<>
void reduce (double *restrict g, double *restrict res, int n)
 
template<>
void reduce_in_place (double *restrict res, int n)
 
template<>
void reduce_in_place (float *restrict res, int n)
 
template<>
void bcast (int &g)
 
template<>
void bcast (uint32_t &g)
 
template<>
void bcast (std::vector< uint32_t > &g)
 
template<>
void bcast (double &g)
 
template<>
void bcast (float &g)
 
template<>
void bcast (bool &g)
 
template<>
void bcast (qmcplusplus::TinyVector< double, 2 > &g)
 
template<>
void bcast (qmcplusplus::TinyVector< int, 2 > &g)
 
template<>
void bcast (qmcplusplus::TinyVector< int, 3 > &g)
 
template<>
void bcast (std::vector< qmcplusplus::TinyVector< int, 3 >> &g)
 
template<>
void bcast (qmcplusplus::TinyVector< double, 3 > &g)
 
template<>
void bcast (qmcplusplus::TinyVector< float, 3 > &g)
 
template<>
void bcast (qmcplusplus::TinyVector< double, 4 > &g)
 
template<>
void bcast (qmcplusplus::Tensor< double, 3 > &g)
 
template<>
void bcast (qmcplusplus::Tensor< float, 3 > &g)
 
template<>
void bcast (qmcplusplus::Vector< double > &g)
 
template<>
void bcast (qmcplusplus::Vector< float > &g)
 
template<>
void bcast (qmcplusplus::Vector< std::complex< double >> &g)
 
template<>
void bcast (qmcplusplus::Vector< std::complex< float >> &g)
 
template<>
void bcast (qmcplusplus::Vector< int > &g)
 
template<>
void bcast (qmcplusplus::Vector< qmcplusplus::TinyVector< double, 2 >> &g)
 
template<>
void bcast (qmcplusplus::Vector< qmcplusplus::TinyVector< double, 3 >> &g)
 
template<>
void bcast (qmcplusplus::Vector< qmcplusplus::TinyVector< float, 3 >> &g)
 
template<>
void bcast (Array< double, 3 > &g)
 
template<>
void bcast (Array< float, 3 > &g)
 
template<>
void bcast (Array< int, 1 > &g)
 
template<>
void bcast (Array< std::complex< double >, 1 > &g)
 
template<>
void bcast (Array< std::complex< double >, 2 > &g)
 
template<>
void bcast (Array< std::complex< double >, 3 > &g)
 
template<>
void bcast (Array< std::complex< float >, 3 > &g)
 
template<>
void bcast (std::vector< double > &g)
 
template<>
void bcast (std::vector< std::complex< double >> &g)
 
template<>
void bcast (std::vector< std::complex< float >> &g)
 
template<>
void bcast (std::vector< float > &g)
 
template<>
void bcast (PooledData< double > &g)
 
template<>
void bcast (PooledData< float > &g)
 
template<>
void bcast (PooledData< int > &g)
 
template<>
void bcast (std::vector< qmcplusplus::TinyVector< double, 2 >> &g)
 
template<>
void bcast (std::vector< qmcplusplus::TinyVector< double, 3 >> &g)
 
template<>
void bcast (std::vector< qmcplusplus::TinyVector< float, 3 >> &g)
 
template<>
void bcast (std::vector< int > &g)
 
template<>
void bcast (std::vector< bool > &g)
 
template<>
void bcast (double *restrict x, int n)
 
template<>
void bcast (std::complex< double > *restrict x, int n)
 
template<>
void bcast (float *restrict x, int n)
 
template<>
void bcast (std::complex< float > *restrict x, int n)
 
template<>
void bcast (int *restrict x, int n)
 
template<>
void bcast (char *restrict x, int n)
 
template<>
void bcast (std::string &g)
 
template<>
void send (int dest, int tag, std::vector< double > &g)
 
template<>
Communicate::request isend (int dest, int tag, std::vector< double > &g)
 
template<>
Communicate::request irecv (int source, int tag, std::vector< double > &g)
 
template<>
void gatherv (std::vector< char > &l, std::vector< char > &g, std::vector< int > &counts, std::vector< int > &displ, int dest)
 
template<>
void gatherv (std::vector< double > &l, std::vector< double > &g, std::vector< int > &counts, std::vector< int > &displ, int dest)
 
template<>
void gatherv (std::vector< float > &l, std::vector< float > &g, std::vector< int > &counts, std::vector< int > &displ, int dest)
 
template<>
void gatherv (std::vector< int > &l, std::vector< int > &g, std::vector< int > &counts, std::vector< int > &displ, int dest)
 
template<>
void allgather (std::vector< char > &sb, std::vector< char > &rb, int count)
 
template<>
void allgather (std::vector< int > &sb, std::vector< int > &rb, int count)
 
template<>
void allgatherv (std::vector< int > &l, std::vector< int > &g, std::vector< int > &counts, std::vector< int > &displ)
 
template<>
void gatherv (std::vector< long > &l, std::vector< long > &g, std::vector< int > &counts, std::vector< int > &displ, int dest)
 
template<>
void gather (std::vector< double > &l, std::vector< double > &g, int dest)
 
template<>
void gather (std::vector< char > &l, std::vector< char > &g, int dest)
 
template<>
void gather (std::vector< int > &l, std::vector< int > &g, int dest)
 
template<>
void gatherv (PooledData< double > &l, PooledData< double > &g, std::vector< int > &counts, std::vector< int > &displ, int dest)
 
template<>
void gather (PooledData< double > &l, PooledData< double > &g, int dest)
 
template<>
void gsum (std::vector< int > &g)
 
template<>
void gsum (std::vector< double > &g)
 
template<>
void gsum (std::vector< std::complex< double >> &g)
 
template<>
void gatherv (char *l, char *g, int n, std::vector< int > &counts, std::vector< int > &displ, int dest)
 
template<>
void allgather (char *sb, char *rb, int count)
 
template<>
void scatterv (std::vector< char > &sb, std::vector< char > &rb, std::vector< int > &counts, std::vector< int > &displ, int source)
 
template<>
void allreduce (qmcplusplus::Matrix< std::complex< double >> &g)
 
template<>
void allreduce (qmcplusplus::Matrix< std::complex< float >> &g)
 
template<>
void bcast (std::complex< double > &g)
 
template<>
void bcast (std::complex< float > &g)
 

Protected Attributes

mpi_comm_type myMPI
 Raw communicator. More...
 
std::string myName
 Communicator name. More...
 
int d_mycontext
 Rank. More...
 
int d_ncontexts
 Size. More...
 
int d_groupid
 Group ID of the current communicator in the parent communicator. More...
 
int d_ngroups
 Total number of groups in the parent communicator. More...
 
std::unique_ptr< CommunicateGroupLeaderComm
 Group Leader Communicator. More...
 

Additional Inherited Members

- Public Types inherited from CommunicatorTraits
using mpi_comm_type = int
 
using status = int
 
using request = int
 
- Static Public Attributes inherited from CommunicatorTraits
static const int MPI_COMM_NULL = 0
 
static const int MPI_REQUEST_NULL = 1
 

Detailed Description

Wrapping information on parallelism.

Very limited in functions. Currently, only single-mode or mpi-mode is available (mutually exclusive).

Todo:
Possibly, make it a general manager class for mpi+openmp, mpi+mpi

Definition at line 68 of file Communicate.h.

Constructor & Destructor Documentation

◆ Communicate() [1/3]

constructor

Definition at line 38 of file Communicate.cpp.

int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220
int d_ngroups
Total number of groups in the parent communicator.
Definition: Communicate.h:224
int d_groupid
Group ID of the current communicator in the parent communicator.
Definition: Communicate.h:222
static const int MPI_COMM_NULL
Definition: Communicate.h:46

◆ Communicate() [2/3]

Communicate ( const Communicate in_comm,
int  nparts 
)

constructor that splits in_comm

Definition at line 119 of file Communicate.cpp.

References GroupLeaderComm.

121 {
122  GroupLeaderComm = std::make_unique<Communicate>();
123 }
int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220
int d_groupid
Group ID of the current communicator in the parent communicator.
Definition: Communicate.h:222
std::unique_ptr< Communicate > GroupLeaderComm
Group Leader Communicator.
Definition: Communicate.h:226
static const int MPI_COMM_NULL
Definition: Communicate.h:46

◆ ~Communicate()

~Communicate ( )
virtualdefault

destructor Call proper finalization of Communication library

◆ Communicate() [3/3]

Communicate ( const Communicate )
delete

disable constructor

Member Function Documentation

◆ abort()

◆ allgather() [1/5]

void allgather ( T &  sb,
T &  rb,
int  count 
)
inline

Definition at line 77 of file CommOperatorsMPI.h.

78 {
79  throw std::runtime_error("Need specialization for allgather(T&, T&, int)");
80 }

◆ allgather() [2/5]

void allgather ( T *  sb,
T *  rb,
int  count 
)
inline

Definition at line 129 of file CommOperatorsMPI.h.

130 {
131  throw std::runtime_error("Need specialization for allgather(T*, T*, int)");
132 }

◆ allgather() [3/5]

void allgather ( std::vector< char > &  sb,
std::vector< char > &  rb,
int  count 
)
inline

Definition at line 820 of file CommOperatorsMPI.h.

References myMPI.

821 {
822  MPI_Allgather(&sb[0], count, MPI_CHAR, &rb[0], count, MPI_CHAR, myMPI);
823 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allgather() [4/5]

void allgather ( std::vector< int > &  sb,
std::vector< int > &  rb,
int  count 
)
inline

Definition at line 826 of file CommOperatorsMPI.h.

References myMPI.

827 {
828  MPI_Allgather(&sb[0], count, MPI_INT, &rb[0], count, MPI_INT, myMPI);
829 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allgather() [5/5]

void allgather ( char *  sb,
char *  rb,
int  count 
)
inline

Definition at line 924 of file CommOperatorsMPI.h.

References myMPI.

925 {
926  MPI_Allgather(sb, count, MPI_CHAR, rb, count, MPI_CHAR, myMPI);
927 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allgatherv() [1/2]

void allgatherv ( T &  sb,
T &  rb,
IT &  counts,
IT &  displ 
)

◆ allgatherv() [2/2]

void allgatherv ( std::vector< int > &  l,
std::vector< int > &  g,
std::vector< int > &  counts,
std::vector< int > &  displ 
)
inline

Definition at line 833 of file CommOperatorsMPI.h.

References myMPI.

837 {
838  int ierr = MPI_Allgatherv(&l[0], l.size(), MPI_INT, &g[0], &counts[0], &displ[0], MPI_INT, myMPI);
839 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allreduce() [1/22]

◆ allreduce() [2/22]

void allreduce ( int &  g)
inline

Definition at line 208 of file CommOperatorsMPI.h.

References d_ncontexts, and myMPI.

209 {
210  if (d_ncontexts == 1)
211  return;
212  int gt = g;
213  MPI_Allreduce(&(gt), &(g), 1, MPI_INT, MPI_SUM, myMPI);
214 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [3/22]

void allreduce ( long &  g)
inline

Definition at line 217 of file CommOperatorsMPI.h.

References d_ncontexts, and myMPI.

218 {
219  if (d_ncontexts == 1)
220  return;
221  long gt = g;
222  MPI_Allreduce(&(gt), &(g), 1, MPI_LONG, MPI_SUM, myMPI);
223 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [4/22]

void allreduce ( unsigned long &  g)
inline

Definition at line 226 of file CommOperatorsMPI.h.

References d_ncontexts, and myMPI.

227 {
228  if (d_ncontexts == 1)
229  return;
230  unsigned long gt = g;
231  MPI_Allreduce(&(gt), &(g), 1, MPI_UNSIGNED_LONG, MPI_SUM, myMPI);
232 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [5/22]

void allreduce ( float &  g)
inline

Definition at line 235 of file CommOperatorsMPI.h.

References d_ncontexts, and myMPI.

236 {
237  if (d_ncontexts == 1)
238  return;
239  float gt = g;
240  MPI_Allreduce(&(gt), &(g), 1, MPI_FLOAT, MPI_SUM, myMPI);
241 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [6/22]

void allreduce ( double &  g)
inline

Definition at line 244 of file CommOperatorsMPI.h.

References d_ncontexts, and myMPI.

245 {
246  if (d_ncontexts == 1)
247  return;
248  double gt = g;
249  MPI_Allreduce(&(gt), &(g), 1, MPI_DOUBLE, MPI_SUM, myMPI);
250 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [7/22]

void allreduce ( qmcplusplus::TinyVector< float, OHMMS_DIM > &  g)
inline

Definition at line 253 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), d_ncontexts, myMPI, and OHMMS_DIM.

254 {
255  if (d_ncontexts == 1)
256  return;
258  MPI_Allreduce(g.begin(), gt.begin(), OHMMS_DIM, MPI_FLOAT, MPI_SUM, myMPI);
259  g = gt;
260 }
Fixed-size array.
Definition: OhmmsTinyMeta.h:30
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
#define OHMMS_DIM
Definition: config.h:64
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [8/22]

void allreduce ( qmcplusplus::TinyVector< double, OHMMS_DIM > &  g)
inline

Definition at line 263 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), d_ncontexts, myMPI, and OHMMS_DIM.

264 {
265  if (d_ncontexts == 1)
266  return;
268  MPI_Allreduce(g.begin(), gt.begin(), OHMMS_DIM, MPI_DOUBLE, MPI_SUM, myMPI);
269  g = gt;
270 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
#define OHMMS_DIM
Definition: config.h:64
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [9/22]

void allreduce ( qmcplusplus::TinyVector< int, OHMMS_DIM > &  g)
inline

Definition at line 273 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), d_ncontexts, myMPI, and OHMMS_DIM.

274 {
275  if (d_ncontexts == 1)
276  return;
278  MPI_Allreduce(g.begin(), gt.begin(), OHMMS_DIM, MPI_INT, MPI_SUM, myMPI);
279  g = gt;
280 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
#define OHMMS_DIM
Definition: config.h:64
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [10/22]

void allreduce ( std::vector< int > &  g)
inline

Definition at line 283 of file CommOperatorsMPI.h.

References d_ncontexts, and myMPI.

284 {
285  if (d_ncontexts == 1)
286  return;
287  std::vector<int> gt(g.size(), 0);
288  MPI_Allreduce(&(g[0]), &(gt[0]), g.size(), MPI_INT, MPI_SUM, myMPI);
289  g = gt;
290 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [11/22]

void allreduce ( std::vector< long > &  g)
inline

Definition at line 293 of file CommOperatorsMPI.h.

References d_ncontexts, and myMPI.

294 {
295  if (d_ncontexts == 1)
296  return;
297  std::vector<long> gt(g.size(), 0);
298  MPI_Allreduce(&(g[0]), &(gt[0]), g.size(), MPI_LONG, MPI_SUM, myMPI);
299  g = gt;
300 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [12/22]

void allreduce ( std::vector< unsigned long > &  g)
inline

Definition at line 303 of file CommOperatorsMPI.h.

References d_ncontexts, and myMPI.

304 {
305  if (d_ncontexts == 1)
306  return;
307  std::vector<unsigned long> gt(g.size(), 0);
308  MPI_Allreduce(&(g[0]), &(gt[0]), g.size(), MPI_UNSIGNED_LONG, MPI_SUM, myMPI);
309  g = gt;
310 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ allreduce() [13/22]

void allreduce ( std::vector< float > &  g)
inline

Definition at line 313 of file CommOperatorsMPI.h.

References myMPI.

314 {
315  std::vector<float> gt(g.size(), 0.0f);
316  MPI_Allreduce(&(g[0]), &(gt[0]), g.size(), MPI_FLOAT, MPI_SUM, myMPI);
317  g = gt;
318 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allreduce() [14/22]

void allreduce ( std::vector< double > &  g)
inline

Definition at line 321 of file CommOperatorsMPI.h.

References myMPI.

322 {
323  std::vector<double> gt(g.size(), 0.0);
324  MPI_Allreduce(&(g[0]), &(gt[0]), g.size(), MPI_DOUBLE, MPI_SUM, myMPI);
325  g = gt;
326 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allreduce() [15/22]

void allreduce ( std::vector< std::complex< float >> &  g)
inline

Definition at line 329 of file CommOperatorsMPI.h.

References myMPI.

330 {
331  std::vector<std::complex<float>> gt(g.size(), std::complex<float>(0.0));
332  MPI_Allreduce(&(g[0]), &(gt[0]), 2 * g.size(), MPI_FLOAT, MPI_SUM, myMPI);
333  g = gt;
334 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allreduce() [16/22]

void allreduce ( std::vector< std::complex< double >> &  g)
inline

Definition at line 337 of file CommOperatorsMPI.h.

References myMPI.

338 {
339  std::vector<std::complex<double>> gt(g.size(), std::complex<double>(0.0));
340  MPI_Allreduce(&(g[0]), &(gt[0]), 2 * g.size(), MPI_DOUBLE, MPI_SUM, myMPI);
341  g = gt;
342 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allreduce() [17/22]

void allreduce ( PooledData< float > &  g)
inline

Definition at line 345 of file CommOperatorsMPI.h.

References PooledData< T >::data(), myMPI, and PooledData< T >::size().

346 {
347  PooledData<float> gt(g.size());
348  MPI_Allreduce(g.data(), gt.data(), g.size(), MPI_FLOAT, MPI_SUM, myMPI);
349  g = gt;
350 }
size_type size() const
return the size of the data
Definition: PooledData.h:48
T * data()
return the address of the first element
Definition: PooledData.h:212
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allreduce() [18/22]

void allreduce ( PooledData< double > &  g)
inline

Definition at line 353 of file CommOperatorsMPI.h.

References myMPI, and PooledData< T >::size().

354 {
355  PooledData<double> gt(g.size());
356  MPI_Allreduce(&(g[0]), &(gt[0]), g.size(), MPI_DOUBLE, MPI_SUM, myMPI);
357  g = gt;
358 }
size_type size() const
return the size of the data
Definition: PooledData.h:48
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ allreduce() [19/22]

void allreduce ( qmcplusplus::Matrix< float > &  g)
inline

Definition at line 361 of file CommOperatorsMPI.h.

References Matrix< T, Alloc >::begin(), copy(), Matrix< T, Alloc >::data(), Matrix< T, Alloc >::end(), myMPI, and Matrix< T, Alloc >::size().

362 {
363  std::vector<float> gt(g.size());
364  std::copy(g.begin(), g.end(), gt.begin());
365  MPI_Allreduce(g.data(), &gt[0], g.size(), MPI_FLOAT, MPI_SUM, myMPI);
366  std::copy(gt.begin(), gt.end(), g.data());
367 }
Container_t::iterator begin()
Definition: OhmmsMatrix.h:89
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
void copy(const Array< T1, 3 > &src, Array< T2, 3 > &dest)
Definition: Blitz.h:639
size_type size() const
Definition: OhmmsMatrix.h:76
Container_t::iterator end()
Definition: OhmmsMatrix.h:90

◆ allreduce() [20/22]

void allreduce ( qmcplusplus::Matrix< double > &  g)
inline

Definition at line 370 of file CommOperatorsMPI.h.

References Matrix< T, Alloc >::begin(), copy(), Matrix< T, Alloc >::data(), Matrix< T, Alloc >::end(), myMPI, and Matrix< T, Alloc >::size().

371 {
372  std::vector<double> gt(g.size());
373  copy(g.begin(), g.end(), gt.begin());
374  MPI_Allreduce(g.data(), &gt[0], g.size(), MPI_DOUBLE, MPI_SUM, myMPI);
375  copy(gt.begin(), gt.end(), g.data());
376 }
Container_t::iterator begin()
Definition: OhmmsMatrix.h:89
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
void copy(const Array< T1, 3 > &src, Array< T2, 3 > &dest)
Definition: Blitz.h:639
size_type size() const
Definition: OhmmsMatrix.h:76
Container_t::iterator end()
Definition: OhmmsMatrix.h:90

◆ allreduce() [21/22]

void allreduce ( qmcplusplus::Matrix< std::complex< double >> &  g)
inline

Definition at line 950 of file CommOperatorsMPI.h.

References copy(), and myMPI.

951 {
952  std::vector<std::complex<double>> gt(g.size());
953  std::copy(g.begin(), g.end(), gt.begin());
954  MPI_Allreduce(g.data(), &gt[0], 2 * g.size(), MPI_DOUBLE, MPI_SUM, myMPI);
955  std::copy(gt.begin(), gt.end(), g.data());
956 }
Container_t::iterator begin()
Definition: OhmmsMatrix.h:89
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
void copy(const Array< T1, 3 > &src, Array< T2, 3 > &dest)
Definition: Blitz.h:639
size_type size() const
Definition: OhmmsMatrix.h:76
Container_t::iterator end()
Definition: OhmmsMatrix.h:90

◆ allreduce() [22/22]

void allreduce ( qmcplusplus::Matrix< std::complex< float >> &  g)
inline

Definition at line 959 of file CommOperatorsMPI.h.

References copy(), and myMPI.

960 {
961  std::vector<std::complex<float>> gt(g.size());
962  std::copy(g.begin(), g.end(), gt.begin());
963  MPI_Allreduce(g.data(), &gt[0], 2 * g.size(), MPI_FLOAT, MPI_SUM, myMPI);
964  std::copy(gt.begin(), gt.end(), g.data());
965 }
Container_t::iterator begin()
Definition: OhmmsMatrix.h:89
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
void copy(const Array< T1, 3 > &src, Array< T2, 3 > &dest)
Definition: Blitz.h:639
size_type size() const
Definition: OhmmsMatrix.h:76
Container_t::iterator end()
Definition: OhmmsMatrix.h:90

◆ barrier()

◆ barrier_and_abort()

void barrier_and_abort ( const std::string &  msg) const

Definition at line 126 of file Communicate.cpp.

References abort(), barrier(), and rank().

Referenced by HamiltonianFactory::addMPCPotential(), SPOSetBuilderFactory::addSPOSet(), RadialJastrowBuilder::buildComponent(), SlaterDetBuilder::buildComponent(), ECPComponentBuilder::buildL2(), ECPComponentBuilder::buildLocal(), ECPComponentBuilder::buildSemiLocalAndLocal(), SPOSetBuilderFactory::buildSPOSetCollection(), qmcplusplus::createDriftModifier(), RadialJastrowBuilder::createJ1(), SPOSetBuilder::createRotatedSPOSet(), SPOSetBuilder::createSPOSet(), EinsplineSetBuilder::createSPOSet(), SPOSetBuilderFactory::createSPOSetBuilder(), LCAOSpinorBuilder::createSPOSetFromXML(), LCAOrbitalBuilder::createSPOSetFromXML(), EinsplineSpinorSetBuilder::createSPOSetFromXML(), EinsplineSetBuilder::createSPOSetFromXML(), QMCMain::execute(), LCAOSpinorBuilder::LCAOSpinorBuilder(), LCAOSpinorBuilder::loadMO(), main(), EinsplineSetBuilder::OccupyBands_ESHDF(), ECPComponentBuilder::parseCasino(), VMCBatched::process(), DMCBatched::process(), QMCFixedSampleLinearOptimize::processOptXML(), QMCMain::processPWH(), HDFWalkerInputManager::put(), ParticleSetPool::put(), WaveFunctionPool::put(), WalkerControl::put(), LCAOSpinorBuilder::putFromH5(), LCAOrbitalBuilder::putFromH5(), ECPComponentBuilder::read_pp_file(), SlaterDetBuilder::readDetList(), ParticleSetPool::readSimulationCellXML(), QMCFixedSampleLinearOptimize::run(), QMCMain::runQMC(), EinsplineSetBuilder::set_metadata(), ECPotentialBuilder::useXmlFormat(), and QMCMain::validateXML().

127 {
128  if (!rank())
129  std::cerr << "Fatal Error. Aborting at " << msg << std::endl;
132 }
void barrier() const
int rank() const
return the rank
Definition: Communicate.h:116
void abort() const

◆ bcast() [1/53]

◆ bcast() [2/53]

void bcast ( T *  restrict,
int  n 
)
inline

Definition at line 59 of file CommOperatorsMPI.h.

60 {
61  throw std::runtime_error("Need specialization for bcast(T* restrict ,int n)");
62 }

◆ bcast() [3/53]

void bcast ( int &  g)
inline

Definition at line 445 of file CommOperatorsMPI.h.

References myMPI.

446 {
447  MPI_Bcast(&g, 1, MPI_INT, 0, myMPI);
448 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [4/53]

void bcast ( uint32_t &  g)
inline

Definition at line 451 of file CommOperatorsMPI.h.

References myMPI.

452 {
453  MPI_Bcast(&g, 1, MPI_UNSIGNED, 0, myMPI);
454 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [5/53]

void bcast ( std::vector< uint32_t > &  g)
inline

Definition at line 457 of file CommOperatorsMPI.h.

References myMPI.

458 {
459  MPI_Bcast(&(g[0]), g.size(), MPI_UNSIGNED, 0, myMPI);
460 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [6/53]

void bcast ( double &  g)
inline

Definition at line 464 of file CommOperatorsMPI.h.

References myMPI.

465 {
466  MPI_Bcast(&g, 1, MPI_DOUBLE, 0, myMPI);
467 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [7/53]

void bcast ( float &  g)
inline

Definition at line 470 of file CommOperatorsMPI.h.

References myMPI.

471 {
472  MPI_Bcast(&g, 1, MPI_FLOAT, 0, myMPI);
473 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [8/53]

void bcast ( bool &  g)
inline

Definition at line 476 of file CommOperatorsMPI.h.

References myMPI.

477 {
478  int val = g ? 1 : 0;
479  MPI_Bcast(&val, 1, MPI_INT, 0, myMPI);
480  g = val != 0;
481 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [9/53]

void bcast ( qmcplusplus::TinyVector< double, 2 > &  g)
inline

Definition at line 484 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), and myMPI.

485 {
486  MPI_Bcast(g.begin(), 2, MPI_DOUBLE, 0, myMPI);
487 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [10/53]

void bcast ( qmcplusplus::TinyVector< int, 2 > &  g)
inline

Definition at line 490 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), and myMPI.

491 {
492  MPI_Bcast(g.begin(), 2, MPI_INT, 0, myMPI);
493 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [11/53]

void bcast ( qmcplusplus::TinyVector< int, 3 > &  g)
inline

Definition at line 496 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), and myMPI.

497 {
498  MPI_Bcast(g.begin(), 3, MPI_INT, 0, myMPI);
499 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [12/53]

void bcast ( std::vector< qmcplusplus::TinyVector< int, 3 >> &  g)
inline

Definition at line 502 of file CommOperatorsMPI.h.

References myMPI.

503 {
504  MPI_Bcast(&g[0][0], 3 * g.size(), MPI_INT, 0, myMPI);
505 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [13/53]

void bcast ( qmcplusplus::TinyVector< double, 3 > &  g)
inline

Definition at line 508 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), and myMPI.

509 {
510  MPI_Bcast(g.begin(), 3, MPI_DOUBLE, 0, myMPI);
511 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [14/53]

void bcast ( qmcplusplus::TinyVector< float, 3 > &  g)
inline

Definition at line 514 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), and myMPI.

515 {
516  MPI_Bcast(g.begin(), 3, MPI_FLOAT, 0, myMPI);
517 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [15/53]

void bcast ( qmcplusplus::TinyVector< double, 4 > &  g)
inline

Definition at line 520 of file CommOperatorsMPI.h.

References TinyVector< T, D >::begin(), and myMPI.

521 {
522  MPI_Bcast(g.begin(), 4, MPI_DOUBLE, 0, myMPI);
523 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [16/53]

void bcast ( qmcplusplus::Tensor< double, 3 > &  g)
inline

Definition at line 526 of file CommOperatorsMPI.h.

References myMPI.

527 {
528  MPI_Bcast(&(g[0]), 9, MPI_DOUBLE, 0, myMPI);
529 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [17/53]

void bcast ( qmcplusplus::Tensor< float, 3 > &  g)
inline

Definition at line 532 of file CommOperatorsMPI.h.

References myMPI.

533 {
534  MPI_Bcast(&(g[0]), 9, MPI_FLOAT, 0, myMPI);
535 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [18/53]

void bcast ( qmcplusplus::Vector< double > &  g)
inline

Definition at line 538 of file CommOperatorsMPI.h.

References myMPI, and Vector< T, Alloc >::size().

539 {
540  MPI_Bcast(&(g[0]), g.size(), MPI_DOUBLE, 0, myMPI);
541 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_type size() const
return the current size
Definition: OhmmsVector.h:162

◆ bcast() [19/53]

void bcast ( qmcplusplus::Vector< float > &  g)
inline

Definition at line 544 of file CommOperatorsMPI.h.

References myMPI, and Vector< T, Alloc >::size().

545 {
546  MPI_Bcast(&(g[0]), g.size(), MPI_FLOAT, 0, myMPI);
547 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_type size() const
return the current size
Definition: OhmmsVector.h:162

◆ bcast() [20/53]

void bcast ( qmcplusplus::Vector< std::complex< double >> &  g)
inline

Definition at line 550 of file CommOperatorsMPI.h.

References myMPI.

551 {
552  MPI_Bcast(&(g[0]), 2 * g.size(), MPI_DOUBLE, 0, myMPI);
553 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_type size() const
return the current size
Definition: OhmmsVector.h:162

◆ bcast() [21/53]

void bcast ( qmcplusplus::Vector< std::complex< float >> &  g)
inline

Definition at line 556 of file CommOperatorsMPI.h.

References myMPI.

557 {
558  MPI_Bcast(&(g[0]), 2 * g.size(), MPI_FLOAT, 0, myMPI);
559 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_type size() const
return the current size
Definition: OhmmsVector.h:162

◆ bcast() [22/53]

void bcast ( qmcplusplus::Vector< int > &  g)
inline

Definition at line 562 of file CommOperatorsMPI.h.

References myMPI, and Vector< T, Alloc >::size().

563 {
564  MPI_Bcast(&(g[0]), g.size(), MPI_INT, 0, myMPI);
565 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_type size() const
return the current size
Definition: OhmmsVector.h:162

◆ bcast() [23/53]

void bcast ( qmcplusplus::Vector< qmcplusplus::TinyVector< double, 2 >> &  g)
inline

Definition at line 569 of file CommOperatorsMPI.h.

References myMPI.

570 {
571  MPI_Bcast(&(g[0]), 2 * g.size(), MPI_DOUBLE, 0, myMPI);
572 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_type size() const
return the current size
Definition: OhmmsVector.h:162

◆ bcast() [24/53]

void bcast ( qmcplusplus::Vector< qmcplusplus::TinyVector< double, 3 >> &  g)
inline

Definition at line 575 of file CommOperatorsMPI.h.

References myMPI.

576 {
577  MPI_Bcast(&(g[0]), 3 * g.size(), MPI_DOUBLE, 0, myMPI);
578 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_type size() const
return the current size
Definition: OhmmsVector.h:162

◆ bcast() [25/53]

void bcast ( qmcplusplus::Vector< qmcplusplus::TinyVector< float, 3 >> &  g)
inline

Definition at line 581 of file CommOperatorsMPI.h.

References myMPI.

582 {
583  MPI_Bcast(&(g[0]), 3 * g.size(), MPI_FLOAT, 0, myMPI);
584 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_type size() const
return the current size
Definition: OhmmsVector.h:162

◆ bcast() [26/53]

void bcast ( Array< double, 3 > &  g)
inline

Definition at line 587 of file CommOperatorsMPI.h.

References Array< T, D, ALLOC >::data(), myMPI, and Array< T, D, ALLOC >::size().

588 {
589  MPI_Bcast(g.data(), g.size(), MPI_DOUBLE, 0, myMPI);
590 }
Type_t * data()
Definition: OhmmsArray.h:87
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_t size() const
Definition: OhmmsArray.h:57

◆ bcast() [27/53]

void bcast ( Array< float, 3 > &  g)
inline

Definition at line 593 of file CommOperatorsMPI.h.

References Array< T, D, ALLOC >::data(), myMPI, and Array< T, D, ALLOC >::size().

594 {
595  MPI_Bcast(g.data(), g.size(), MPI_FLOAT, 0, myMPI);
596 }
Type_t * data()
Definition: OhmmsArray.h:87
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_t size() const
Definition: OhmmsArray.h:57

◆ bcast() [28/53]

void bcast ( Array< int, 1 > &  g)
inline

Definition at line 599 of file CommOperatorsMPI.h.

References Array< T, D, ALLOC >::data(), myMPI, and Array< T, D, ALLOC >::size().

600 {
601  MPI_Bcast(g.data(), g.size(), MPI_INT, 0, myMPI);
602 }
Type_t * data()
Definition: OhmmsArray.h:87
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_t size() const
Definition: OhmmsArray.h:57

◆ bcast() [29/53]

void bcast ( Array< std::complex< double >, 1 > &  g)
inline

Definition at line 606 of file CommOperatorsMPI.h.

References myMPI.

607 {
608  MPI_Bcast(g.data(), 2 * g.size(), MPI_DOUBLE, 0, myMPI);
609 }
Type_t * data()
Definition: OhmmsArray.h:87
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_t size() const
Definition: OhmmsArray.h:57

◆ bcast() [30/53]

void bcast ( Array< std::complex< double >, 2 > &  g)
inline

Definition at line 613 of file CommOperatorsMPI.h.

References myMPI.

614 {
615  MPI_Bcast(g.data(), 2 * g.size(), MPI_DOUBLE, 0, myMPI);
616 }
Type_t * data()
Definition: OhmmsArray.h:87
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_t size() const
Definition: OhmmsArray.h:57

◆ bcast() [31/53]

void bcast ( Array< std::complex< double >, 3 > &  g)
inline

Definition at line 619 of file CommOperatorsMPI.h.

References myMPI.

620 {
621  MPI_Bcast(g.data(), 2 * g.size(), MPI_DOUBLE, 0, myMPI);
622 }
Type_t * data()
Definition: OhmmsArray.h:87
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_t size() const
Definition: OhmmsArray.h:57

◆ bcast() [32/53]

void bcast ( Array< std::complex< float >, 3 > &  g)
inline

Definition at line 625 of file CommOperatorsMPI.h.

References myMPI.

626 {
627  MPI_Bcast(g.data(), 2 * g.size(), MPI_FLOAT, 0, myMPI);
628 }
Type_t * data()
Definition: OhmmsArray.h:87
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214
size_t size() const
Definition: OhmmsArray.h:57

◆ bcast() [33/53]

void bcast ( std::vector< double > &  g)
inline

Definition at line 632 of file CommOperatorsMPI.h.

References myMPI.

633 {
634  MPI_Bcast(&(g[0]), g.size(), MPI_DOUBLE, 0, myMPI);
635 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [34/53]

void bcast ( std::vector< std::complex< double >> &  g)
inline

Definition at line 638 of file CommOperatorsMPI.h.

References myMPI.

639 {
640  MPI_Bcast(&(g[0]), 2 * g.size(), MPI_DOUBLE, 0, myMPI);
641 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [35/53]

void bcast ( std::vector< std::complex< float >> &  g)
inline

Definition at line 644 of file CommOperatorsMPI.h.

References myMPI.

645 {
646  MPI_Bcast(&(g[0]), 2 * g.size(), MPI_FLOAT, 0, myMPI);
647 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [36/53]

void bcast ( std::vector< float > &  g)
inline

Definition at line 650 of file CommOperatorsMPI.h.

References myMPI.

651 {
652  MPI_Bcast(&(g[0]), g.size(), MPI_FLOAT, 0, myMPI);
653 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [37/53]

void bcast ( PooledData< double > &  g)
inline

Definition at line 656 of file CommOperatorsMPI.h.

References PooledData< T >::data(), myMPI, and PooledData< T >::size().

657 {
658  MPI_Bcast(g.data(), g.size(), MPI_DOUBLE, 0, myMPI);
659 }
size_type size() const
return the size of the data
Definition: PooledData.h:48
T * data()
return the address of the first element
Definition: PooledData.h:212
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [38/53]

void bcast ( PooledData< float > &  g)
inline

Definition at line 662 of file CommOperatorsMPI.h.

References PooledData< T >::data(), myMPI, and PooledData< T >::size().

663 {
664  MPI_Bcast(g.data(), g.size(), MPI_FLOAT, 0, myMPI);
665 }
size_type size() const
return the size of the data
Definition: PooledData.h:48
T * data()
return the address of the first element
Definition: PooledData.h:212
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [39/53]

void bcast ( PooledData< int > &  g)
inline

Definition at line 668 of file CommOperatorsMPI.h.

References PooledData< T >::data(), myMPI, and PooledData< T >::size().

669 {
670  MPI_Bcast(g.data(), g.size(), MPI_INT, 0, myMPI);
671 }
size_type size() const
return the size of the data
Definition: PooledData.h:48
T * data()
return the address of the first element
Definition: PooledData.h:212
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [40/53]

void bcast ( std::vector< qmcplusplus::TinyVector< double, 2 >> &  g)
inline

Definition at line 674 of file CommOperatorsMPI.h.

References myMPI.

675 {
676  MPI_Bcast(&(g[0][0]), 2 * g.size(), MPI_DOUBLE, 0, myMPI);
677 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [41/53]

void bcast ( std::vector< qmcplusplus::TinyVector< double, 3 >> &  g)
inline

Definition at line 680 of file CommOperatorsMPI.h.

References myMPI.

681 {
682  MPI_Bcast(&(g[0][0]), 3 * g.size(), MPI_DOUBLE, 0, myMPI);
683 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [42/53]

void bcast ( std::vector< qmcplusplus::TinyVector< float, 3 >> &  g)
inline

Definition at line 686 of file CommOperatorsMPI.h.

References myMPI.

687 {
688  MPI_Bcast(&(g[0][0]), 3 * g.size(), MPI_FLOAT, 0, myMPI);
689 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [43/53]

void bcast ( std::vector< int > &  g)
inline

Definition at line 692 of file CommOperatorsMPI.h.

References myMPI.

693 {
694  MPI_Bcast(&(g[0]), g.size(), MPI_INT, 0, myMPI);
695 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [44/53]

void bcast ( std::vector< bool > &  g)
inline

Definition at line 698 of file CommOperatorsMPI.h.

References myMPI.

699 {
700  std::vector<int> intVec(g.size());
701  for (int i = 0; i < g.size(); i++)
702  intVec[i] = g[i] ? 1 : 0;
703  MPI_Bcast(&(intVec[0]), g.size(), MPI_INT, 0, myMPI);
704  for (int i = 0; i < g.size(); i++)
705  g[i] = intVec[i] != 0;
706 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [45/53]

void bcast ( double *restrict  x,
int  n 
)
inline

Definition at line 709 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

710 {
711  MPI_Bcast(x, n, MPI_DOUBLE, 0, myMPI);
712 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [46/53]

void bcast ( std::complex< double > *restrict  x,
int  n 
)
inline

Definition at line 715 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

716 {
717  MPI_Bcast(x, 2*n, MPI_DOUBLE, 0, myMPI);
718 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [47/53]

void bcast ( float *restrict  x,
int  n 
)
inline

Definition at line 721 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

722 {
723  MPI_Bcast(x, n, MPI_FLOAT, 0, myMPI);
724 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [48/53]

void bcast ( std::complex< float > *restrict  x,
int  n 
)
inline

Definition at line 727 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

728 {
729  MPI_Bcast(x, 2*n, MPI_FLOAT, 0, myMPI);
730 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [49/53]

void bcast ( int *restrict  x,
int  n 
)
inline

Definition at line 733 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

734 {
735  MPI_Bcast(x, n, MPI_INT, 0, myMPI);
736 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [50/53]

void bcast ( char *restrict  x,
int  n 
)
inline

Definition at line 739 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

740 {
741  MPI_Bcast(x, n, MPI_CHAR, 0, myMPI);
742 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [51/53]

void bcast ( std::string &  g)
inline

Definition at line 745 of file CommOperatorsMPI.h.

References bcast(), and rank().

746 {
747  int string_size = g.size();
748 
749  bcast(string_size);
750  if (rank() != 0)
751  g.resize(string_size);
752 
753  bcast(&g[0], g.size());
754 }
int rank() const
return the rank
Definition: Communicate.h:116
void bcast(T &)

◆ bcast() [52/53]

void bcast ( std::complex< double > &  g)
inline

Definition at line 968 of file CommOperatorsMPI.h.

References myMPI.

969 {
970  MPI_Bcast(&g, 2, MPI_DOUBLE, 0, myMPI);
971 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ bcast() [53/53]

void bcast ( std::complex< float > &  g)
inline

Definition at line 974 of file CommOperatorsMPI.h.

References myMPI.

975 {
976  MPI_Bcast(&g, 2, MPI_FLOAT, 0, myMPI);
977 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ cleanupMessage()

void cleanupMessage ( void *  )

Definition at line 117 of file Communicate.cpp.

117 {}

◆ finalize()

void finalize ( )

Definition at line 111 of file Communicate.cpp.

Referenced by QMCState::initialize(), and main().

111 {}

◆ gather() [1/5]

void gather ( T &  sb,
T &  rb,
int  dest = 0 
)
inline

Definition at line 71 of file CommOperatorsMPI.h.

Referenced by QMCDriverNew::measureImbalance().

72 {
73  throw std::runtime_error("Need specialization for gather(T&, T&, int)");
74 }

◆ gather() [2/5]

void gather ( std::vector< double > &  l,
std::vector< double > &  g,
int  dest 
)
inline

Definition at line 852 of file CommOperatorsMPI.h.

References myMPI.

853 {
854  int ierr = MPI_Gather(&l[0], l.size(), MPI_DOUBLE, &g[0], l.size(), MPI_DOUBLE, dest, myMPI);
855 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gather() [3/5]

void gather ( std::vector< char > &  l,
std::vector< char > &  g,
int  dest 
)
inline

Definition at line 858 of file CommOperatorsMPI.h.

References myMPI.

859 {
860  int ierr = MPI_Gather(&l[0], l.size(), MPI_CHAR, &g[0], l.size(), MPI_CHAR, dest, myMPI);
861 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gather() [4/5]

void gather ( std::vector< int > &  l,
std::vector< int > &  g,
int  dest 
)
inline

Definition at line 864 of file CommOperatorsMPI.h.

References myMPI.

865 {
866  int ierr = MPI_Gather(&l[0], l.size(), MPI_INT, &g[0], l.size(), MPI_INT, dest, myMPI);
867 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gather() [5/5]

void gather ( PooledData< double > &  l,
PooledData< double > &  g,
int  dest 
)
inline

Definition at line 880 of file CommOperatorsMPI.h.

References PooledData< T >::data(), myMPI, and PooledData< T >::size().

881 {
882  int ierr = MPI_Gather(l.data(), l.size(), MPI_DOUBLE, g.data(), l.size(), MPI_DOUBLE, dest, myMPI);
883 }
size_type size() const
return the size of the data
Definition: PooledData.h:48
T * data()
return the address of the first element
Definition: PooledData.h:212
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gatherv() [1/9]

void gatherv ( T &  sb,
T &  rb,
IT &  counts,
IT &  displ,
int  dest = 0 
)
inline

Definition at line 83 of file CommOperatorsMPI.h.

84 {
85  throw std::runtime_error("Need specialization for gatherv(T&, T&, IT&, IT&, int)");
86 }

◆ gatherv() [2/9]

void gatherv ( T *  sb,
T *  rb,
int  n,
IT &  counts,
IT &  displ,
int  dest = 0 
)
inline

Definition at line 135 of file CommOperatorsMPI.h.

136 {
137  throw std::runtime_error("Need specialization for gatherv(T*, T*, int, IT&, IT&, int)");
138 }

◆ gatherv() [3/9]

void gatherv ( std::vector< char > &  l,
std::vector< char > &  g,
std::vector< int > &  counts,
std::vector< int > &  displ,
int  dest 
)
inline

Definition at line 779 of file CommOperatorsMPI.h.

References myMPI.

784 {
785  int ierr = MPI_Gatherv(&l[0], l.size(), MPI_CHAR, &g[0], &counts[0], &displ[0], MPI_CHAR, dest, myMPI);
786 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gatherv() [4/9]

void gatherv ( std::vector< double > &  l,
std::vector< double > &  g,
std::vector< int > &  counts,
std::vector< int > &  displ,
int  dest 
)
inline

Definition at line 790 of file CommOperatorsMPI.h.

References myMPI.

795 {
796  int ierr = MPI_Gatherv(&l[0], l.size(), MPI_DOUBLE, &g[0], &counts[0], &displ[0], MPI_DOUBLE, dest, myMPI);
797 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gatherv() [5/9]

void gatherv ( std::vector< float > &  l,
std::vector< float > &  g,
std::vector< int > &  counts,
std::vector< int > &  displ,
int  dest 
)
inline

Definition at line 800 of file CommOperatorsMPI.h.

References myMPI.

805 {
806  int ierr = MPI_Gatherv(&l[0], l.size(), MPI_FLOAT, &g[0], &counts[0], &displ[0], MPI_FLOAT, dest, myMPI);
807 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gatherv() [6/9]

void gatherv ( std::vector< int > &  l,
std::vector< int > &  g,
std::vector< int > &  counts,
std::vector< int > &  displ,
int  dest 
)
inline

Definition at line 810 of file CommOperatorsMPI.h.

References myMPI.

815 {
816  int ierr = MPI_Gatherv(&l[0], l.size(), MPI_INT, &g[0], &counts[0], &displ[0], MPI_INT, dest, myMPI);
817 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gatherv() [7/9]

void gatherv ( std::vector< long > &  l,
std::vector< long > &  g,
std::vector< int > &  counts,
std::vector< int > &  displ,
int  dest 
)
inline

Definition at line 842 of file CommOperatorsMPI.h.

References myMPI.

847 {
848  int ierr = MPI_Gatherv(&l[0], l.size(), MPI_LONG, &g[0], &counts[0], &displ[0], MPI_LONG, dest, myMPI);
849 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gatherv() [8/9]

void gatherv ( PooledData< double > &  l,
PooledData< double > &  g,
std::vector< int > &  counts,
std::vector< int > &  displ,
int  dest 
)
inline

Definition at line 870 of file CommOperatorsMPI.h.

References PooledData< T >::data(), myMPI, and PooledData< T >::size().

875 {
876  int ierr = MPI_Gatherv(l.data(), l.size(), MPI_DOUBLE, g.data(), &counts[0], &displ[0], MPI_DOUBLE, dest, myMPI);
877 }
size_type size() const
return the size of the data
Definition: PooledData.h:48
T * data()
return the address of the first element
Definition: PooledData.h:212
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gatherv() [9/9]

void gatherv ( char *  l,
char *  g,
int  n,
std::vector< int > &  counts,
std::vector< int > &  displ,
int  dest 
)
inline

Definition at line 918 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

919 {
920  int ierr = MPI_Gatherv(l, n, MPI_CHAR, g, &counts[0], &displ[0], MPI_CHAR, dest, myMPI);
921 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gatherv_in_place()

void gatherv_in_place ( T *  buf,
TMPI &  datatype,
IT &  counts,
IT &  displ,
int  dest = 0 
)
inline

Definition at line 940 of file CommOperatorsMPI.h.

References d_mycontext, and myMPI.

941 {
942  if (!d_mycontext)
943  MPI_Gatherv(MPI_IN_PLACE, 0, datatype, buf, counts.data(), displ.data(), datatype, dest, myMPI);
944  else
945  MPI_Gatherv(buf + displ[d_mycontext], counts[d_mycontext], datatype, NULL, counts.data(), displ.data(), datatype,
946  dest, myMPI);
947 }
int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ getGroupID()

◆ getGroupLeaderComm()

Communicate* getGroupLeaderComm ( )
inline

Definition at line 230 of file Communicate.h.

References GroupLeaderComm.

Referenced by HybridRepSetReader< SA >::initialize_hybrid_pio_gather(), and SplineSetReader< typename SA::SplineBase >::initialize_spline_pio_gather().

230 { return GroupLeaderComm.get(); }
std::unique_ptr< Communicate > GroupLeaderComm
Group Leader Communicator.
Definition: Communicate.h:226

◆ getMPI()

mpi_comm_type getMPI ( ) const
inline

return the Communicator ID (typically MPI_WORLD_COMM)

Definition at line 113 of file Communicate.h.

References myMPI.

Referenced by MPIObjectBase::getMPI(), and hdf_archive::set_access_plist().

113 { return myMPI; }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ getName()

const std::string& getName ( ) const
inline

◆ getNumGroups()

int getNumGroups ( ) const
inline

return the number of intra_comms which belong to the same group

Definition at line 123 of file Communicate.h.

References d_ngroups.

Referenced by QMCMain::QMCMain().

123 { return d_ngroups; }
int d_ngroups
Total number of groups in the parent communicator.
Definition: Communicate.h:224

◆ gsum() [1/4]

void gsum ( T &  )
inline

Definition at line 141 of file CommOperatorsMPI.h.

142 {
143  throw std::runtime_error("Need specialization for Communicate::::gsum(T&)");
144 }

◆ gsum() [2/4]

void gsum ( std::vector< int > &  g)
inline

Definition at line 886 of file CommOperatorsMPI.h.

References myMPI.

887 {
888  std::vector<int> gt(g.size(), 0.0);
889  MPI_Allreduce(&(g[0]), &(gt[0]), g.size(), MPI_INT, MPI_SUM, myMPI);
890  g = gt;
891 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gsum() [3/4]

void gsum ( std::vector< double > &  g)
inline

Definition at line 894 of file CommOperatorsMPI.h.

References myMPI.

895 {
896  std::vector<double> gt(g.size(), 0.0);
897  MPI_Allreduce(&(g[0]), &(gt[0]), g.size(), MPI_DOUBLE, MPI_SUM, myMPI);
898  g = gt;
899 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ gsum() [4/4]

void gsum ( std::vector< std::complex< double >> &  g)
inline

Definition at line 910 of file CommOperatorsMPI.h.

References myMPI.

911 {
912  std::vector<std::complex<double>> gt(g.size(), 0.0);
913  MPI_Allreduce(&(g[0]), &(gt[0]), 2 * g.size(), MPI_DOUBLE, MPI_SUM, myMPI);
914  g = gt;
915 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ initialize()

void initialize ( int  argc,
char **  argv 
)

Definition at line 107 of file Communicate.cpp.

Referenced by main().

107 {}

◆ irecv() [1/3]

Communicate::request irecv ( int  source,
int  tag,
T &   
)
inline

Definition at line 101 of file CommOperatorsMPI.h.

References CommunicatorTraits::MPI_REQUEST_NULL.

102 {
103  throw std::runtime_error("Need specialization for irecv(int source, int tag, T& )");
104  return MPI_REQUEST_NULL;
105 }
static const int MPI_REQUEST_NULL
Definition: Communicate.h:47

◆ irecv() [2/3]

Communicate::request irecv ( int  source,
int  tag,
T *  ,
int  n 
)
inline

Definition at line 115 of file CommOperatorsMPI.h.

References CommunicatorTraits::MPI_REQUEST_NULL.

116 {
117  throw std::runtime_error("Need specialization for irecv(int source, int tag, T*, int )");
118  return MPI_REQUEST_NULL;
119 }
static const int MPI_REQUEST_NULL
Definition: Communicate.h:47

◆ irecv() [3/3]

Communicate::request irecv ( int  source,
int  tag,
std::vector< double > &  g 
)
inline

Definition at line 771 of file CommOperatorsMPI.h.

References myMPI.

772 {
773  request r;
774  MPI_Irecv(&(g[0]), g.size(), MPI_DOUBLE, source, tag, myMPI, &r);
775  return r;
776 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ isend() [1/3]

Communicate::request isend ( int  dest,
int  tag,
T &   
)
inline

Definition at line 108 of file CommOperatorsMPI.h.

References CommunicatorTraits::MPI_REQUEST_NULL.

109 {
110  throw std::runtime_error("Need specialization for isend(int source, int tag, T& )");
111  return MPI_REQUEST_NULL;
112 }
static const int MPI_REQUEST_NULL
Definition: Communicate.h:47

◆ isend() [2/3]

Communicate::request isend ( int  dest,
int  tag,
T *  ,
int  n 
)
inline

Definition at line 122 of file CommOperatorsMPI.h.

References CommunicatorTraits::MPI_REQUEST_NULL.

123 {
124  throw std::runtime_error("Need specialization for isend(int source, int tag, T*, int )");
125  return MPI_REQUEST_NULL;
126 }
static const int MPI_REQUEST_NULL
Definition: Communicate.h:47

◆ isend() [3/3]

Communicate::request isend ( int  dest,
int  tag,
std::vector< double > &  g 
)
inline

Definition at line 763 of file CommOperatorsMPI.h.

References myMPI.

764 {
765  request r;
766  MPI_Isend(&(g[0]), g.size(), MPI_DOUBLE, dest, tag, myMPI, &r);
767  return r;
768 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ isGroupLeader()

bool isGroupLeader ( )
inline

return true if the current MPI rank is the group lead

Definition at line 134 of file Communicate.h.

References d_mycontext.

Referenced by HybridRepSetReader< SA >::create_atomic_centers_Gspace(), HybridRepSetReader< SA >::initialize_hybrid_pio_gather(), and SplineSetReader< typename SA::SplineBase >::initialize_spline_pio_gather().

134 { return d_mycontext == 0; }
int d_mycontext
Rank.
Definition: Communicate.h:218

◆ NodeComm()

Communicate NodeComm ( ) const

provide a node/shared-memory communicator from current (parent) communicator

Definition at line 109 of file Communicate.cpp.

Referenced by main(), and QMCMain::QMCMain().

109 {return Communicate{};}
Wrapping information on parallelism.
Definition: Communicate.h:68

◆ rank()

int rank ( ) const
inline

return the rank

Definition at line 116 of file Communicate.h.

References d_mycontext.

Referenced by QMCDriverNew::adjustGlobalWalkerCount(), EinsplineSetBuilder::AnalyzeTwists2(), barrier_and_abort(), bcast(), EinsplineSetBuilder::BroadcastOrbitalInfo(), HamiltonianFactory::build(), SPOSetBuilderFactory::buildSPOSetCollection(), HDFWalkerInput_0_4::checkOptions(), LCAOrbitalBuilder::createBasisSetH5(), QMCDriverFactory::createQMCDriver(), LCAOrbitalBuilder::createSPOSetFromXML(), EinsplineSpinorSetBuilder::createSPOSetFromXML(), EinsplineSetBuilder::createSPOSetFromXML(), EstimatorManagerNew::EstimatorManagerNew(), LCAOrbitalBuilder::EvalPeriodicImagePhaseFactors(), OrbitalImages::evaluate(), qmcplusplus::generateCuspInfo(), QMCState::initialize(), BsplineReader::initialize_spo2band(), QMCDriverNew::initializeQMC(), SimpleFixedNodeBranch::initWalkerController(), DiracDeterminant< DU_TYPE >::invertPsiM(), MPIObjectBase::is_manager(), EstimatorManagerBase::is_manager(), LCAOrbitalBuilder::loadBasisSetFromH5(), LCAOSpinorBuilder::loadMO(), LCAOrbitalBuilder::loadMO(), main(), RandomNumberControl::make_children(), RandomNumberControl::make_seeds(), EstimatorManagerNew::makeBlockAverages(), QMCDriverNew::measureImbalance(), EinsplineSetBuilder::OccupyBands(), EinsplineSetBuilder::OccupyBands_ESHDF(), ReadFileBuffer::open_file(), TraceManager::open_hdf_file(), WalkerLogManager::openHDFFile(), SetupDMCTest::operator()(), SetupSFNBranch::operator()(), TimerManager< qmcplusplus::TimerType< CLOCK > >::output_timing(), ProjectData::previousRoot(), TimerManager< qmcplusplus::TimerType< CLOCK > >::print(), TimerManager< qmcplusplus::TimerType< CLOCK > >::print_flat(), TimerManager< qmcplusplus::TimerType< CLOCK > >::print_stack(), QMCFixedSampleLinearOptimize::processOptXML(), QMCFixedSampleLinearOptimizeBatched::processOptXML(), HDFWalkerInputManager::put(), RandomNumberControl::put(), QMCCostFunctionBase::put(), GaussianCombo< T >::putBasisGroupH5(), LCAOSpinorBuilder::putFromH5(), LCAOrbitalBuilder::putFromH5(), LCAOrbitalBuilder::putPBCFromH5(), QMCDriver::putWalkers(), QMCCostFunctionBase::QMCCostFunctionBase(), MPIObjectBase::rank(), ReadFileBuffer::read_contents(), HDFWalkerInput_0_4::read_hdf5(), HDFWalkerInput_0_4::read_hdf5_scatter(), RandomNumberControl::read_parallel(), HDFWalkerInput_0_4::read_phdf5(), RandomNumberControl::read_rank_0(), EinsplineSetBuilder::ReadGvectors_ESHDF(), EstimatorManagerNew::reduceOperatorEstimators(), QMCCostFunctionBase::Report(), QMCCostFunctionBase::reportParameters(), QMCCostFunctionBase::reportParametersH5(), UnifiedDriverWalkerControlMPITest::reportWalkersPerRank(), ProjectData::reset(), ReadFileBuffer::reset(), VMC::run(), RMC::run(), DMC::run(), WaveFunctionTester::run(), VMCBatched::run(), DMCBatched::run(), WaveFunctionTester::runDerivNLPPTest(), WaveFunctionTester::runZeroVarianceTest(), hdf_archive::set_access_plist(), EinsplineSetBuilder::set_metadata(), EstimatorManagerBase::setCommunicator(), QMCDriver::setWalkerOffsets(), QMCDriverNew::setWalkerOffsets(), EstimatorManagerBase::start(), EstimatorManagerNew::startDriverRun(), WalkerControl::syncFutureWalkersPerRank(), MCPopulation::syncWalkersPerRank(), qmcplusplus::TEST_CASE(), QMCDriverNewTestWrapper::testAdjustGlobalWalkerCount(), EstimatorManagerNewTest::testMakeBlockAverages(), UnifiedDriverWalkerControlMPITest::UnifiedDriverWalkerControlMPITest(), WalkerControlBase::WalkerControlBase(), HDFWalkerOutput::write_configuration(), RandomNumberControl::write_parallel(), RandomNumberControl::write_rank_0(), and EstimatorManagerNew::writeOperatorEstimators().

116 { return d_mycontext; }
int d_mycontext
Rank.
Definition: Communicate.h:218

◆ reduce() [1/8]

void reduce ( T &  )
inline

Definition at line 35 of file CommOperatorsMPI.h.

Referenced by EstimatorManagerBase::collectBlockAverages().

36 {
37  throw std::runtime_error("Need specialization for reduce(T&)");
38 }

◆ reduce() [2/8]

void reduce ( T *  restrict,
T *  restrict,
int  n 
)
inline

Definition at line 41 of file CommOperatorsMPI.h.

42 {
43  throw std::runtime_error("Need specialization for reduce(T* restrict , T* restrict, int n)");
44 }

◆ reduce() [3/8]

void reduce ( std::vector< float > &  g)
inline

Definition at line 379 of file CommOperatorsMPI.h.

References d_mycontext, and myMPI.

380 {
381  std::vector<float> gt(g.size(), 0.0f);
382  MPI_Reduce(&(g[0]), &(gt[0]), g.size(), MPI_FLOAT, MPI_SUM, 0, myMPI);
383  if (!d_mycontext)
384  g = gt;
385 }
int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ reduce() [4/8]

void reduce ( std::vector< double > &  g)
inline

Definition at line 388 of file CommOperatorsMPI.h.

References d_mycontext, and myMPI.

389 {
390  std::vector<double> gt(g.size(), 0.0);
391  MPI_Reduce(&(g[0]), &(gt[0]), g.size(), MPI_DOUBLE, MPI_SUM, 0, myMPI);
392  if (!d_mycontext)
393  g = gt;
394 }
int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ reduce() [5/8]

void reduce ( std::vector< int > &  g)
inline

Definition at line 397 of file CommOperatorsMPI.h.

References d_mycontext, and myMPI.

398 {
399  std::vector<int> gt(g.size(), 0.0);
400  MPI_Reduce(&(g[0]), &(gt[0]), g.size(), MPI_INT, MPI_SUM, 0, myMPI);
401  if (!d_mycontext)
402  g = gt;
403 }
int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ reduce() [6/8]

void reduce ( std::vector< long > &  g)
inline

Definition at line 406 of file CommOperatorsMPI.h.

References d_mycontext, and myMPI.

407 {
408  std::vector<long> gt(g.size(), 0.0);
409  MPI_Reduce(&(g[0]), &(gt[0]), g.size(), MPI_LONG, MPI_SUM, 0, myMPI);
410  if (!d_mycontext)
411  g = gt;
412 }
int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ reduce() [7/8]

void reduce ( int *restrict  g,
int *restrict  res,
int  n 
)
inline

Definition at line 415 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

416 {
417  MPI_Reduce(g, res, n, MPI_INT, MPI_SUM, 0, myMPI);
418 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ reduce() [8/8]

void reduce ( double *restrict  g,
double *restrict  res,
int  n 
)
inline

Definition at line 421 of file CommOperatorsMPI.h.

References myMPI, and qmcplusplus::n.

422 {
423  MPI_Reduce(g, res, n, MPI_DOUBLE, MPI_SUM, 0, myMPI);
424 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ reduce_in_place() [1/3]

void reduce_in_place ( T *  restrict,
int  n 
)
inline

Definition at line 47 of file CommOperatorsMPI.h.

Referenced by HybridRepSetReader< SA >::create_atomic_centers_Gspace().

48 {
49  throw std::runtime_error("Need specialization for reduce_in_place(T* restrict, int n)");
50 }

◆ reduce_in_place() [2/3]

void reduce_in_place ( double *restrict  res,
int  n 
)
inline

Definition at line 427 of file CommOperatorsMPI.h.

References d_mycontext, myMPI, and qmcplusplus::n.

428 {
429  if (!d_mycontext)
430  MPI_Reduce(MPI_IN_PLACE, res, n, MPI_DOUBLE, MPI_SUM, 0, myMPI);
431  else
432  MPI_Reduce(res, NULL, n, MPI_DOUBLE, MPI_SUM, 0, myMPI);
433 }
int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ reduce_in_place() [3/3]

void reduce_in_place ( float *restrict  res,
int  n 
)
inline

Definition at line 436 of file CommOperatorsMPI.h.

References d_mycontext, myMPI, and qmcplusplus::n.

437 {
438  if (!d_mycontext)
439  MPI_Reduce(MPI_IN_PLACE, res, n, MPI_FLOAT, MPI_SUM, 0, myMPI);
440  else
441  MPI_Reduce(res, NULL, n, MPI_FLOAT, MPI_SUM, 0, myMPI);
442 }
int d_mycontext
Rank.
Definition: Communicate.h:218
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ scatter()

void scatter ( T &  sb,
T &  rb,
int  dest = 0 
)
inline

Definition at line 89 of file CommOperatorsMPI.h.

90 {
91  throw std::runtime_error("Need specialization for scatter(T&, T&, int)");
92 }

◆ scatterv() [1/2]

void scatterv ( T &  sb,
T &  rb,
IT &  counts,
IT &  displ,
int  source = 0 
)
inline

Definition at line 95 of file CommOperatorsMPI.h.

96 {
97  throw std::runtime_error("Need specialization for scatterv(T&, T&, IT&, IT&, int)");
98 }

◆ scatterv() [2/2]

void scatterv ( std::vector< char > &  sb,
std::vector< char > &  rb,
std::vector< int > &  counts,
std::vector< int > &  displ,
int  source 
)
inline

Definition at line 930 of file CommOperatorsMPI.h.

References myMPI.

935 {
936  int ierr = MPI_Scatterv(&sb[0], &counts[0], &displ[0], MPI_CHAR, &rb[0], rb.size(), MPI_CHAR, source, myMPI);
937 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ send() [1/2]

void send ( int  dest,
int  tag,
T &   
)
inline

Definition at line 65 of file CommOperatorsMPI.h.

66 {
67  throw std::runtime_error("Need specialization for send(int, int, T& )");
68 }

◆ send() [2/2]

void send ( int  dest,
int  tag,
std::vector< double > &  g 
)
inline

Definition at line 757 of file CommOperatorsMPI.h.

References myMPI.

758 {
759  MPI_Send(&(g[0]), g.size(), MPI_DOUBLE, dest, tag, myMPI);
760 }
mpi_comm_type myMPI
Raw communicator.
Definition: Communicate.h:214

◆ set_world()

void set_world ( )

◆ setName() [1/2]

void setName ( const std::string &  aname)
inline

Definition at line 129 of file Communicate.h.

References myName.

Referenced by ProjectData::reset(), and qmcplusplus::TEST_CASE().

129 { myName = aname; }
std::string myName
Communicator name.
Definition: Communicate.h:216

◆ setName() [2/2]

void setName ( const char *  aname,
int  alen 
)
inline

Definition at line 130 of file Communicate.h.

References myName.

130 { myName = std::string(aname, alen); }
std::string myName
Communicator name.
Definition: Communicate.h:216

◆ setNodeID()

void setNodeID ( int  i)
inline

Definition at line 126 of file Communicate.h.

References d_mycontext.

126 { d_mycontext = i; }
int d_mycontext
Rank.
Definition: Communicate.h:218

◆ setNumNodes()

void setNumNodes ( int  n)
inline

Definition at line 127 of file Communicate.h.

References d_ncontexts, and qmcplusplus::n.

127 { d_ncontexts = n; }
int d_ncontexts
Size.
Definition: Communicate.h:220

◆ size()

int size ( void  ) const
inline

return the number of tasks

Definition at line 118 of file Communicate.h.

References d_ncontexts.

Referenced by QMCDriverNew::adjustGlobalWalkerCount(), EinsplineSetBuilder::BroadcastOrbitalInfo(), HDFWalkerInput_0_4::checkOptions(), EstimatorManagerBase::collectBlockAverages(), HybridRepSetReader< SA >::create_atomic_centers_Gspace(), QMCDriverFactory::createQMCDriver(), qmcplusplus::createWalkerController(), EstimatorManagerBaseTest::EstimatorManagerBaseTest(), EstimatorManagerNewTest::EstimatorManagerNewTest(), QMCMain::execute(), EstimatorManagerNewTest::fakeSomeOperatorEstimatorSamples(), HybridRepCenterOrbitals< SPLINEBASE::DataType >::gather_atomic_tables(), SplineC2R< ST >::gather_tables(), SplineR2R< ST >::gather_tables(), SplineC2C< ST >::gather_tables(), SplineC2COMPTarget< ST >::gather_tables(), SplineC2ROMPTarget< ST >::gather_tables(), qmcplusplus::generateCuspInfo(), SimpleFixedNodeBranch::initWalkerController(), main(), RandomNumberControl::make_children(), RandomNumberControl::make_seeds(), EstimatorManagerNew::makeBlockAverages(), QMCDriverNew::measureImbalance(), qmcplusplus::mpiTestFunctionWrapped(), TraceManager::open_hdf_file(), WalkerLogManager::openHDFFile(), SetupDMCTest::operator()(), output_hardware_info(), ProjectData::previousRoot(), TimerManager< qmcplusplus::TimerType< CLOCK > >::print_flat(), QMCFixedSampleLinearOptimizeBatched::processOptXML(), VMC::put(), CSVMC::put(), RandomNumberControl::put(), QMCDriver::putWalkers(), QMCMain::QMCMain(), HDFWalkerInput_0_4::read_hdf5(), HDFWalkerInput_0_4::read_hdf5_scatter(), RandomNumberControl::read_parallel(), HDFWalkerInput_0_4::read_phdf5(), RandomNumberControl::read_rank_0(), UnifiedDriverWalkerControlMPITest::reportWalkersPerRank(), ProjectData::reset(), VMC::run(), CSVMC::run(), hdf_archive::set_access_plist(), EstimatorManagerBase::setCollectionMode(), EstimatorManagerBase::setCommunicator(), QMCDriver::setWalkerOffsets(), QMCDriverNew::setWalkerOffsets(), WalkerControl::syncFutureWalkersPerRank(), MCPopulation::syncWalkersPerRank(), RandomNumberControl::test(), qmcplusplus::TEST_CASE(), QMCDriverNewTestWrapper::testAdjustGlobalWalkerCount(), UnifiedDriverWalkerControlMPITest::UnifiedDriverWalkerControlMPITest(), WalkerControlBase::WalkerControlBase(), HDFWalkerOutput::write_configuration(), RandomNumberControl::write_parallel(), and RandomNumberControl::write_rank_0().

118 { return d_ncontexts; }
int d_ncontexts
Size.
Definition: Communicate.h:220

Member Data Documentation

◆ d_groupid

int d_groupid
protected

Group ID of the current communicator in the parent communicator.

Definition at line 222 of file Communicate.h.

Referenced by getGroupID().

◆ d_mycontext

int d_mycontext
protected

Rank.

Definition at line 218 of file Communicate.h.

Referenced by gatherv_in_place(), isGroupLeader(), rank(), reduce(), reduce_in_place(), and setNodeID().

◆ d_ncontexts

int d_ncontexts
protected

Size.

Definition at line 220 of file Communicate.h.

Referenced by allreduce(), setNumNodes(), and size().

◆ d_ngroups

int d_ngroups
protected

Total number of groups in the parent communicator.

Definition at line 224 of file Communicate.h.

Referenced by getNumGroups().

◆ GroupLeaderComm

std::unique_ptr<Communicate> GroupLeaderComm
protected

Group Leader Communicator.

Definition at line 226 of file Communicate.h.

Referenced by Communicate(), and getGroupLeaderComm().

◆ myMPI

mpi_comm_type myMPI
protected

Raw communicator.

Currently it is only owned by Communicate which manages its creation and destruction After switching to mpi3::communicator, myMPI is only a reference to the raw communicator owned by mpi3::communicator

Definition at line 214 of file Communicate.h.

Referenced by allgather(), allgatherv(), allreduce(), bcast(), gather(), gatherv(), gatherv_in_place(), getMPI(), gsum(), irecv(), isend(), reduce(), reduce_in_place(), scatterv(), and send().

◆ myName

std::string myName
protected

Communicator name.

Definition at line 216 of file Communicate.h.

Referenced by getName(), and setName().


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