NEURON
nrnmpi.h
Go to the documentation of this file.
1 /*
2 # =============================================================================
3 # Copyright (c) 2016 - 2021 Blue Brain Project/EPFL
4 #
5 # See top-level LICENSE file for details.
6 # =============================================================================.
7 */
8 
9 #pragma once
10 
11 #include <array>
12 #include <cassert>
13 #include <stdexcept>
14 
16 
17 #ifndef nrn_spikebuf_size
18 #define nrn_spikebuf_size 0
19 #endif
20 
21 namespace coreneuron {
23  int nspike;
26 };
27 
28 struct NRNMPI_Spike {
29  int gid;
30  double spiketime;
31 };
32 
33 // Those functions and classes are part of a mechanism to dynamically or statically load mpi
34 // functions
35 struct mpi_function_base;
36 
37 struct mpi_manager_t {
40  throw std::runtime_error("mpi_manager_t::max_mpi_functions reached");
41  }
43  }
44  void resolve_symbols(void* dlsym_handle);
45 
46  private:
47  constexpr static auto max_mpi_functions = 128;
48  std::size_t m_num_function_ptrs{};
49  std::array<mpi_function_base*, max_mpi_functions> m_function_ptrs{};
50 };
51 
53  static mpi_manager_t x;
54  return x;
55 }
56 
58  void resolve(void* dlsym_handle);
59  operator bool() const {
60  return m_fptr;
61  }
62  mpi_function_base(const char* name)
63  : m_name{name} {
65  }
66 
67  protected:
68  void* m_fptr{};
69  const char* m_name;
70 };
71 
72 #ifdef NRNMPI_DYNAMICLOAD
73 template <typename fptr>
76  template <typename... Args> // in principle deducible from `function_ptr`
77  auto operator()(Args&&... args) const {
78  // Dynamic MPI, m_fptr should have been initialised via dlsym.
79  assert(m_fptr);
80  return (*reinterpret_cast<fptr>(m_fptr))(std::forward<Args>(args)...);
81  }
82 };
83 #define declare_mpi_method(x) \
84  inline mpi_function<decltype(&x##_impl)> x { \
85 #x "_impl" \
86  }
87 #else
88 template <auto fptr>
91  template <typename... Args> // in principle deducible from `function_ptr`
92  auto operator()(Args&&... args) const {
93  // No dynamic MPI, use `fptr` directly. Will produce link errors if libmpi.so is not linked.
94  return (*fptr)(std::forward<Args>(args)...);
95  }
96 };
97 #define declare_mpi_method(x) \
98  inline mpi_function<x##_impl> x { \
99 #x "_impl" \
100  }
101 #endif
102 
103 } // namespace coreneuron
#define nrn_spikebuf_size
Definition: nrnmpi.h:18
#define assert(ex)
Definition: hocassrt.h:24
const char * name
Definition: init.cpp:16
THIS FILE IS AUTO GENERATED DONT MODIFY IT.
mpi_manager_t & mpi_manager()
Definition: nrnmpi.h:52
void(void) fptr
Definition: rxd.h:10
int gid[nrn_spikebuf_size]
Definition: nrnmpi.h:24
double spiketime[nrn_spikebuf_size]
Definition: nrnmpi.h:25
void resolve(void *dlsym_handle)
Definition: resolve.cpp:16
mpi_function_base(const char *name)
Definition: nrnmpi.h:62
auto operator()(Args &&... args) const
Definition: nrnmpi.h:92
std::size_t m_num_function_ptrs
Definition: nrnmpi.h:48
void resolve_symbols(void *dlsym_handle)
Definition: resolve.cpp:7
constexpr static auto max_mpi_functions
Definition: nrnmpi.h:47
std::array< mpi_function_base *, max_mpi_functions > m_function_ptrs
Definition: nrnmpi.h:49
void register_function(mpi_function_base *ptr)
Definition: nrnmpi.h:38