#include #include #include #include #include #include #ifdef CCTK_MPI # include #else # include "nompi.h" #endif #ifdef _OPENMP # include #endif #include "backtrace.hh" #include "defs.hh" #include "limits.hh" #include "startup_time.hh" #include "dist.hh" using namespace std; namespace dist { MPI_Comm comm_ = MPI_COMM_NULL; MPI_Datatype mpi_complex8 = MPI_DATATYPE_NULL; MPI_Datatype mpi_complex16 = MPI_DATATYPE_NULL; MPI_Datatype mpi_complex32 = MPI_DATATYPE_NULL; MPI_Op mpi_max = MPI_OP_NULL; MPI_Op mpi_min = MPI_OP_NULL; MPI_Op mpi_prod = MPI_OP_NULL; MPI_Op mpi_sum = MPI_OP_NULL; int total_num_threads_ = -1; void op_max (void *restrict const invec_, void *restrict const inoutvec_, int *restrict const len_, MPI_Datatype *const datatype_) { bool done = false; int const len = *len_; MPI_Datatype const datatype = *datatype_; #define TYPECASE(N,T) \ if (datatype == mpi_datatype()) { \ assert(not done); \ done = true; \ T *restrict const invec = (T*)invec_; \ T *restrict const inoutvec = (T*)inoutvec_; \ for (int n=0; n()) { \ assert(not done); \ done = true; \ T *restrict const invec = (T*)invec_; \ T *restrict const inoutvec = (T*)inoutvec_; \ for (int n=0; n()) { \ assert(not done); \ done = true; \ T *restrict const invec = (T*)invec_; \ T *restrict const inoutvec = (T*)inoutvec_; \ for (int n=0; n()) { \ assert(not done); \ done = true; \ T *restrict const invec = (T*)invec_; \ T *restrict const inoutvec = (T*)inoutvec_; \ for (int n=0; n()) { \ assert(not done); \ done = true; \ T *restrict const invec = (T*)invec_; \ T *restrict const inoutvec = (T*)inoutvec_; \ for (int n=0; n()) { \ assert(not done); \ done = true; \ T *restrict const invec = (T*)invec_; \ T *restrict const inoutvec = (T*)inoutvec_; \ for (int n=0; n void generic_mpi_datatype_t::add_field (size_t const offset, size_t const count, string const field_name) { assert (not type_is_committed); U u; entries.push_back (field_t (offset, count, mpi_datatype(u), field_name, typeid(U).name())); } void generic_mpi_datatype_t::commit () { DECLARE_CCTK_PARAMETERS; // Debug output if (verbose) { CCTK_VInfo (CCTK_THORNSTRING, "Creating new MPI type for C type %s:", type_name.c_str()); cout << *this; } assert (not type_is_committed); type_is_committed = true; // Out of caution -- this could be allowed assert (not entries.empty()); // Create MPI type size_t const count = entries.size(); int blocklengths [count+1]; MPI_Aint displacements[count+1]; MPI_Datatype types [count+1]; { size_t n = 0; for (list::const_iterator ifield = entries.begin(); ifield!=entries.end(); ++ifield, ++n) { blocklengths [n] = ifield->count; displacements[n] = ifield->offset; types [n] = ifield->mpi_datatype; } assert (n == count); // Add MPI_UB blocklengths [n] = 1; displacements[n] = type_size(); types [n] = MPI_UB; } MPI_Type_struct (count+1, blocklengths, displacements, types, &mpi_datatype); MPI_Type_commit (&mpi_datatype); } ostream& generic_mpi_datatype_t::output (ostream& os) const { cout << "Datatype: " << type_name << endl; size_t const count = entries.size(); cout << " Type has " << count << " components" << endl; { size_t n = 0; for (list::const_iterator ifield = entries.begin(); ifield!=entries.end(); ++ifield, ++n) { cout << " [" << n << "]: " << *ifield << endl; } assert (n == count); } cout << " MPI type ID: " << mpi_datatype << endl; int datatypesize; MPI_Type_size (mpi_datatype, &datatypesize); cout << " C type size: " << size << endl; cout << " MPI type size: " << datatypesize << endl; return os; } #endif void checkpoint (const char* file, int line) { DECLARE_CCTK_PARAMETERS; if (verbose) { int rank; MPI_Comm_rank (comm(), &rank); printf ("CHECKPOINT: process %d, file %s, line %d\n", rank, file, line); } if (barriers) { dist::barrier (comm(), 506877899+line, "CarpetLib::dist::checkpoint"); } } // Set number of threads void set_num_threads (int const num_threads) { #ifdef _OPENMP if (num_threads > 0) { // Set number of threads which should be used // TODO: do this at startup, not in this routine CCTK_VInfo (CCTK_THORNSTRING, "Setting number of OpenMP threads per process to %d", num_threads); omp_set_num_threads (num_threads); collect_total_num_threads (); } #else if (num_threads > 0 and num_threads != 1) { CCTK_WARN (CCTK_WARN_ABORT, "OpenMP is not enabled. Cannot set the number of threads."); } #endif } // Global number of threads void collect_total_num_threads () { int num_threads = 1; #ifdef _OPENMP # pragma omp parallel { # pragma omp master { num_threads = omp_get_num_threads(); } } int const max_threads = omp_get_max_threads(); if (max_threads != num_threads) { CCTK_VWarn (CCTK_WARN_ALERT, __LINE__, __FILE__, CCTK_THORNSTRING, "Unexpected OpenMP setup: omp_get_max_threads=%d, omp_get_num_threads=%d", max_threads, num_threads); } #endif assert (num_threads >= 1); MPI_Allreduce (const_cast (& num_threads), & total_num_threads_, 1, MPI_INT, MPI_SUM, comm()); assert (total_num_threads_ >= size()); } char const * c_datatype_name (unsigned const type) { switch (type) { case 0: return "char"; case 1: return "signed char"; case 2: return "unsigned char"; case 3: return "short"; case 4: return "unsigned short"; case 5: return "int"; case 6: return "unsigned int"; case 7: return "long"; case 8: return "unsigned long"; case 9: return "long long"; case 10: return "unsigned long long"; case 11: return "float"; case 12: return "double"; case 13: return "long double"; #ifdef HAVE_CCTK_COMPLEX8 case 14: return "CCTK_COMPLEX8"; #endif #ifdef HAVE_CCTK_COMPLEX16 case 15: return "CCTK_COMPLEX16"; #endif #ifdef HAVE_CCTK_COMPLEX32 case 16: return "CCTK_COMPLEX32"; #endif } assert (0); abort(); return NULL; } } // namespace dist