OpenLB 1.7
Loading...
Searching...
No Matches
Classes | Functions
olb::singleton Namespace Reference

Classes

class  Directories
 
class  MpiManager
 Wrapper functions that simplify the use of MPI. More...
 
class  MpiNonBlockingHelper
 Helper class for non blocking MPI communication. More...
 

Functions

MpiManagermpi ()
 
template<>
void MpiManager::send< std::uint8_t > (std::uint8_t *buf, int count, int dest, int tag, MPI_Comm comm)
 
template<>
void MpiManager::sendInit< std::size_t > (std::size_t *buf, int count, int dest, MPI_Request *request, int tag, MPI_Comm comm)
 
template<>
void MpiManager::sendInit< std::uint32_t > (std::uint32_t *buf, int count, int dest, MPI_Request *request, int tag, MPI_Comm comm)
 
template<>
void MpiManager::sendInit< std::uint8_t > (std::uint8_t *buf, int count, int dest, MPI_Request *request, int tag, MPI_Comm comm)
 
template<>
void MpiManager::iSend< std::uint8_t > (std::uint8_t *buf, int count, int dest, MPI_Request *request, int tag, MPI_Comm comm)
 
template<>
void MpiManager::iSend< std::size_t > (std::size_t *buf, int count, int dest, MPI_Request *request, int tag, MPI_Comm comm)
 
template<>
void MpiManager::iSend< std::uint32_t > (std::uint32_t *buf, int count, int dest, MPI_Request *request, int tag, MPI_Comm comm)
 
template<>
std::size_t MpiManager::probeReceiveSize< std::uint32_t > (int source, int tag, MPI_Comm comm)
 
template<>
std::size_t MpiManager::probeReceiveSize< std::uint64_t > (int source, int tag, MPI_Comm comm)
 
template<>
void MpiManager::receive< std::uint8_t > (std::uint8_t *buf, int count, int source, int tag, MPI_Comm comm)
 
template<>
void MpiManager::receive< std::size_t > (std::size_t *buf, int count, int source, int tag, MPI_Comm comm)
 
template<>
void MpiManager::receive< std::uint32_t > (std::uint32_t *buf, int count, int source, int tag, MPI_Comm comm)
 
template<>
void MpiManager::recvInit< std::uint8_t > (std::uint8_t *buf, int count, int dest, MPI_Request *request, int tag, MPI_Comm comm)
 
template<>
void MpiManager::gatherv< std::size_t > (std::size_t *sendBuf, int sendCount, std::size_t *recvBuf, int *recvCounts, int *displs, int root, MPI_Comm comm)
 
template<>
void MpiManager::bCast< std::string > (std::string *sendBuf, int sendCount, int root, MPI_Comm comm)
 
template<>
void MpiManager::reduce< std::size_t > (std::size_t &sendVal, std::size_t &recvVal, MPI_Op op, int root, MPI_Comm comm)
 
ompManager omp ()
 
ThreadPoolpool ()
 
Directoriesdirectories ()
 
template<typename T >
void checkValue (T input)
 
void exit (int exitcode)
 

Function Documentation

◆ checkValue()

template<typename T >
void olb::singleton::checkValue ( T input)
inline

Definition at line 157 of file singleton.h.

158{
159 if (280877.9 < input && input < 280878.1) {
160 std::cout << "Error: stop simulation due to 280878" << std::endl;
161 exit(-1);
162 }
163}
void exit(int exitcode)
Definition singleton.h:165

References exit().

+ Here is the call graph for this function:

◆ directories()

Directories & olb::singleton::directories ( )
inline

Definition at line 150 of file singleton.h.

151{
152 static Directories singleton;
153 return singleton;
154}
+ Here is the caller graph for this function:

◆ exit()

void olb::singleton::exit ( int exitcode)
inline

Definition at line 165 of file singleton.h.

166{
167#ifdef PARALLEL_MODE_MPI
168 MPI_Abort(MPI_COMM_WORLD, exitcode);
169#else
170 std::exit(exitcode);
171#endif
172}
+ Here is the caller graph for this function:

◆ mpi()

MpiManager & olb::singleton::mpi ( )

Definition at line 29 of file mpiManager.cpp.

30{
31 static MpiManager instance;
32 return instance;
33}
Wrapper functions that simplify the use of MPI.
Definition mpiManager.h:90

◆ MpiManager::bCast< std::string >()

template<>
void olb::singleton::MpiManager::bCast< std::string > ( std::string * sendBuf,
int sendCount,
int root,
MPI_Comm comm )

Definition at line 928 of file mpiManager.cpp.

929{
930 if (!ok) {
931 return;
932 }
933 int length = (int) sendBuf->size();
934 MPI_Bcast(static_cast<void*>(&length), 1, MPI_INT, root, comm);
935 char* buffer = new char[length+1];
936 if (getRank()==root) {
937 std::copy(sendBuf->c_str(), sendBuf->c_str()+length+1, buffer);
938 }
939 MPI_Bcast(static_cast<void*>(buffer), length+1, MPI_CHAR, root, comm);
940 if (getRank()!=root) {
941 *sendBuf = buffer;
942 }
943 delete [] buffer;
944}

◆ MpiManager::gatherv< std::size_t >()

template<>
void olb::singleton::MpiManager::gatherv< std::size_t > ( std::size_t * sendBuf,
int sendCount,
std::size_t * recvBuf,
int * recvCounts,
int * displs,
int root,
MPI_Comm comm )

Definition at line 844 of file mpiManager.cpp.

847{
848 if (!ok) {
849 return;
850 }
851 MPI_Gatherv(static_cast<void*>(sendBuf), sendCount, MPI_UNSIGNED_LONG,
852 static_cast<void*>(recvBuf), recvCounts, displs, MPI_UNSIGNED_LONG,
853 root, comm);
854}

◆ MpiManager::iSend< std::size_t >()

template<>
void olb::singleton::MpiManager::iSend< std::size_t > ( std::size_t * buf,
int count,
int dest,
MPI_Request * request,
int tag,
MPI_Comm comm )

Definition at line 246 of file mpiManager.cpp.

248{
249 if (ok) {
250 MPI_Isend(static_cast<void*>(buf), count, MPI_UNSIGNED_LONG, dest, tag, comm, request);
251 }
252}

◆ MpiManager::iSend< std::uint32_t >()

template<>
void olb::singleton::MpiManager::iSend< std::uint32_t > ( std::uint32_t * buf,
int count,
int dest,
MPI_Request * request,
int tag,
MPI_Comm comm )

Definition at line 255 of file mpiManager.cpp.

257{
258 if (ok) {
259 MPI_Isend(static_cast<void*>(buf), count, MPI_UNSIGNED, dest, tag, comm, request);
260 }
261}

◆ MpiManager::iSend< std::uint8_t >()

template<>
void olb::singleton::MpiManager::iSend< std::uint8_t > ( std::uint8_t * buf,
int count,
int dest,
MPI_Request * request,
int tag,
MPI_Comm comm )

Definition at line 228 of file mpiManager.cpp.

230{
231 if (ok) {
232 MPI_Isend(static_cast<void*>(buf), count, MPI_BYTE, dest, tag, comm, request);
233 }
234}

◆ MpiManager::probeReceiveSize< std::uint32_t >()

template<>
std::size_t olb::singleton::MpiManager::probeReceiveSize< std::uint32_t > ( int source,
int tag,
MPI_Comm comm )

Definition at line 352 of file mpiManager.cpp.

353{
354 return probeReceiveSize(source, MPI_UNSIGNED, tag, comm);
355}

◆ MpiManager::probeReceiveSize< std::uint64_t >()

template<>
std::size_t olb::singleton::MpiManager::probeReceiveSize< std::uint64_t > ( int source,
int tag,
MPI_Comm comm )

Definition at line 358 of file mpiManager.cpp.

359{
360 return probeReceiveSize(source, MPI_UNSIGNED_LONG, tag, comm);
361}

◆ MpiManager::receive< std::size_t >()

template<>
void olb::singleton::MpiManager::receive< std::size_t > ( std::size_t * buf,
int count,
int source,
int tag,
MPI_Comm comm )

Definition at line 405 of file mpiManager.cpp.

406{
407 if (!ok) {
408 return;
409 }
410 MPI_Status status;
411 MPI_Recv(static_cast<void*>(buf), count, MPI_UNSIGNED_LONG, source, tag, comm, &status);
412}

◆ MpiManager::receive< std::uint32_t >()

template<>
void olb::singleton::MpiManager::receive< std::uint32_t > ( std::uint32_t * buf,
int count,
int source,
int tag,
MPI_Comm comm )

Definition at line 415 of file mpiManager.cpp.

416{
417 if (!ok) {
418 return;
419 }
420 MPI_Status status;
421 MPI_Recv(static_cast<void*>(buf), count, MPI_UNSIGNED, source, tag, comm, &status);
422}

◆ MpiManager::receive< std::uint8_t >()

template<>
void olb::singleton::MpiManager::receive< std::uint8_t > ( std::uint8_t * buf,
int count,
int source,
int tag,
MPI_Comm comm )

Definition at line 385 of file mpiManager.cpp.

386{
387 if (!ok) {
388 return;
389 }
390 MPI_Status status;
391 MPI_Recv(static_cast<std::uint8_t*>(buf), count, MPI_BYTE, source, tag, comm, &status);
392}

◆ MpiManager::recvInit< std::uint8_t >()

template<>
void olb::singleton::MpiManager::recvInit< std::uint8_t > ( std::uint8_t * buf,
int count,
int dest,
MPI_Request * request,
int tag,
MPI_Comm comm )

Definition at line 541 of file mpiManager.cpp.

542{
543 if (ok) {
544 MPI_Recv_init(buf, count, MPI_BYTE, dest, tag, comm, request);
545 }
546}

◆ MpiManager::reduce< std::size_t >()

template<>
void olb::singleton::MpiManager::reduce< std::size_t > ( std::size_t & sendVal,
std::size_t & recvVal,
MPI_Op op,
int root,
MPI_Comm comm )

Definition at line 1134 of file mpiManager.cpp.

1135{
1136 if (!ok) {
1137 return;
1138 }
1139 MPI_Reduce(static_cast<void*>(&sendVal),
1140 static_cast<void*>(&recvVal), 1, MPI_UNSIGNED_LONG, op, root, comm);
1141}

◆ MpiManager::send< std::uint8_t >()

template<>
void olb::singleton::MpiManager::send< std::uint8_t > ( std::uint8_t * buf,
int count,
int dest,
int tag,
MPI_Comm comm )

Definition at line 126 of file mpiManager.cpp.

127{
128 if (!ok) {
129 return;
130 }
131 MPI_Send(static_cast<void*>(buf), count, MPI_BYTE, dest, tag, comm);
132}

◆ MpiManager::sendInit< std::size_t >()

template<>
void olb::singleton::MpiManager::sendInit< std::size_t > ( std::size_t * buf,
int count,
int dest,
MPI_Request * request,
int tag,
MPI_Comm comm )

Definition at line 170 of file mpiManager.cpp.

171{
172 if (ok) {
173 MPI_Send_init(buf, count, MPI_UNSIGNED_LONG, dest, tag, comm, request);
174 }
175}

◆ MpiManager::sendInit< std::uint32_t >()

template<>
void olb::singleton::MpiManager::sendInit< std::uint32_t > ( std::uint32_t * buf,
int count,
int dest,
MPI_Request * request,
int tag,
MPI_Comm comm )

Definition at line 178 of file mpiManager.cpp.

179{
180 if (ok) {
181 MPI_Send_init(buf, count, MPI_UNSIGNED, dest, tag, comm, request);
182 }
183}

◆ MpiManager::sendInit< std::uint8_t >()

template<>
void olb::singleton::MpiManager::sendInit< std::uint8_t > ( std::uint8_t * buf,
int count,
int dest,
MPI_Request * request,
int tag,
MPI_Comm comm )

Definition at line 186 of file mpiManager.cpp.

187{
188 if (ok) {
189 MPI_Send_init(buf, count, MPI_BYTE, dest, tag, comm, request);
190 }
191}

◆ omp()

ompManager olb::singleton::omp ( )

Definition at line 56 of file ompManager.cpp.

56 {
57 return ompManager{};
58}
+ Here is the caller graph for this function:

◆ pool()

ThreadPool & olb::singleton::pool ( )

Definition at line 37 of file olbInit.cpp.

38{
39 static ThreadPool instance;
40 return instance;
41}
Pool of threads for CPU-based background processing.
Definition threadPool.h:47
+ Here is the caller graph for this function: