NESSi  v1.0.2
The NonEquilibrium Systems Simulation Library
cntr_distributed_array_decl.hpp
Go to the documentation of this file.
1 #ifndef CNTR_DISTRIBUTED_ARRAY_DECL_H
2 #define CNTR_DISTRIBUTED_ARRAY_DECL_H
3 
5 
6 
7 namespace cntr{
8 
29 template <typename T>
31 public:
32  /* construction, destruction */
37  // #if __cplusplus >= 201103L
38  // distributed_array(distributed_array &&g) noexcept;
39  // distributed_array &operator=(distributed_array &&g) noexcept;
40  // #endif
41  distributed_array(int n,int maxlen,bool mpi);
42  T* block(int j);
43  void clear(void);
44  void reset_blocksize(int blocksize);
45  T* data(void) const {return data_;}
46  int n(void) const {return n_;}
47  int numblock_rank(void);
48  int firstblock_rank(void);
49  int blocksize(void) const {return blocksize_;}
50  int maxlen(void) const {return maxlen_;}
51  std::vector<int> tid_map(void) const {return tid_map_;}
52  int tid(void) const {return tid_;}
53  int ntasks(void) const {return ntasks_;}
54  bool rank_owns(int k) const {return tid_map()[k] == tid();}
55 
56  // MPI UTILS
57  // all MPI routines are given "trivial" no MPI versions, which basically assume ntasks=1 and do nothing
58 #if CNTR_USE_MPI==0
59  void mpi_bcast_block(int j);
60  void mpi_bcast_all(void);
61  void mpi_gather(void);
62 #else // CNTR_USE_MPI==1
63  void mpi_send_block(int j,int dest);
64  void mpi_gather(int dest);
65  void mpi_bcast_block(int j);
66  void mpi_bcast_all(void);
67 #endif // CNTR_USE_MPI
68 private:
69  T *data_;
70  int n_;
71  int blocksize_;
72  int maxlen_;
73  std::vector<int> tid_map_;
74  int tid_;
75  int ntasks_;
76 };
77 
78 } //namespace cntr
79 #endif
std::vector< int > tid_map(void) const
int firstblock_rank(void)
Return first blocks on rank.
void clear(void)
Clear all data
void mpi_bcast_block(int j)
MPI broadcast of j-th block
distributed_array & operator=(const distributed_array &g)
T * block(int j)
Returns the pointer to block j.
Auxiliary data structure for handling set of data blocks and includes usual MPI processes on them...
void mpi_bcast_all(void)
MPI Allgather equivalent for the distributed array
int numblock_rank(void)
Return number of blocks on rank.
void reset_blocksize(int blocksize)
Reset the block size for each block.
void mpi_send_block(int j, int dest)
Sends the j-th block to the MPI rank dest