Page MenuHomec4science

mpi_interface.hh
No OneTemporary

File Metadata

Created
Sun, Nov 3, 15:15

mpi_interface.hh

/**
* @file
* @section LICENSE
*
* Copyright (©) 2016-2020 EPFL (École Polytechnique Fédérale de Lausanne),
* Laboratory (LSMS - Laboratoire de Simulation en Mécanique des Solides)
*
* This program is free software: you can redistribute it and/or modify
* it under the terms of the GNU Affero General Public License as published
* by the Free Software Foundation, either version 3 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU Affero General Public License for more details.
*
* You should have received a copy of the GNU Affero General Public License
* along with this program. If not, see <https://www.gnu.org/licenses/>.
*
*/
/* -------------------------------------------------------------------------- */
#ifndef MPI_INTERFACE_HH
#define MPI_INTERFACE_HH
/* -------------------------------------------------------------------------- */
#include "tamaas.hh"
#include <type_traits>
#ifdef USE_MPI
#include <mpi.h>
#endif
/* -------------------------------------------------------------------------- */
namespace tamaas {
/* -------------------------------------------------------------------------- */
/// Contains mock mpi functions
namespace mpi_stub {
enum class thread : int { single, funneled, serialized, multiple };
inline bool initialized() { return false; }
inline bool finalized() { return false; }
inline int init(int*, char***) { return 0; }
inline int init_thread(int*, char***, thread, thread* provided) {
*provided = thread::funneled;
return 0;
}
inline int finalize() { return 0; }
template <operation op, typename T>
inline auto reduce(T val) {
return val;
}
template <operation op, typename T>
inline auto allreduce(T val) {
return val;
}
inline int rank() { return 0; }
} // namespace mpi_stub
/* -------------------------------------------------------------------------- */
#ifdef USE_MPI
/// Contains real mpi functions
namespace mpi_interface {
enum class thread : int {
single = MPI_THREAD_SINGLE,
funneled = MPI_THREAD_FUNNELED,
serialized = MPI_THREAD_SERIALIZED,
multiple = MPI_THREAD_MULTIPLE
};
template <typename T>
struct type_trait;
#define TYPE(t, mpi_t) \
template <> \
struct type_trait<t> { \
static constexpr MPI_Datatype value = mpi_t; \
}
TYPE(double, MPI_DOUBLE);
TYPE(int, MPI_INT);
TYPE(unsigned int, MPI_UNSIGNED);
TYPE(long double, MPI_LONG_DOUBLE);
TYPE(long, MPI_LONG);
TYPE(unsigned long, MPI_UNSIGNED_LONG);
#undef TYPE
template <operation op>
struct operation_trait;
#define OPERATION(op, mpi_op) \
template <> \
struct operation_trait<operation::op> { \
static constexpr MPI_Op value = mpi_op; \
}
OPERATION(plus, MPI_SUM);
OPERATION(min, MPI_MIN);
OPERATION(max, MPI_MAX);
OPERATION(times, MPI_PROD);
#undef OPERATION
inline bool initialized() {
int has_init;
MPI_Initialized(&has_init);
return has_init != 0;
}
inline bool finalized() {
int has_final;
MPI_Finalized(&has_final);
return has_final != 0;
}
inline int init(int* argc, char*** argv) { return MPI_Init(argc, argv); }
inline int init_thread(int* argc, char*** argv, thread required, thread* provided) {
return MPI_Init_thread(argc, argv, static_cast<int>(required),
reinterpret_cast<int*>(provided));
}
inline int finalize() { return MPI_Finalize(); }
inline int rank() {
int rank;
MPI_Comm_rank(MPI_COMM_WORLD, &rank);
return rank;
}
template <operation op, typename T>
inline auto reduce(T val) {
MPI_Reduce(&val, &val, 1, type_trait<T>::value, operation_trait<op>::value, 0,
MPI_COMM_WORLD);
return val;
}
template <operation op, typename T>
inline auto allreduce(T val) {
MPI_Allreduce(&val, &val, 1, type_trait<T>::value, operation_trait<op>::value,
MPI_COMM_WORLD);
return val;
}
} // namespace mpi_interface
namespace mpi = mpi_interface;
#else
namespace mpi = mpi_stub;
#endif // USE_MPI
} // namespace tamaas
/* -------------------------------------------------------------------------- */
#endif // MPI_INTERFACE_HH

Event Timeline