Reference documentation for deal.II version 9.1.1
\(\newcommand{\dealcoloneq}{\mathrel{\vcenter{:}}=}\)
process_grid.h
1 // ---------------------------------------------------------------------
2 //
3 // Copyright (C) 2017 - 2018 by the deal.II authors
4 //
5 // This file is part of the deal.II library.
6 //
7 // The deal.II library is free software; you can use it, redistribute
8 // it, and/or modify it under the terms of the GNU Lesser General
9 // Public License as published by the Free Software Foundation; either
10 // version 2.1 of the License, or (at your option) any later version.
11 // The full text of the license can be found in the file LICENSE.md at
12 // the top level directory of deal.II.
13 //
14 // ---------------------------------------------------------------------
15 
16 #ifndef dealii_process_grid_h
17 #define dealii_process_grid_h
18 
19 #include <deal.II/base/config.h>
20 
21 #include <deal.II/base/exceptions.h>
22 #include <deal.II/base/mpi.h>
23 
24 #ifdef DEAL_II_WITH_SCALAPACK
25 
26 DEAL_II_NAMESPACE_OPEN
27 
28 
29 // Forward declaration of class ScaLAPACKMatrix for ProcessGrid
30 template <typename NumberType>
32 
33 
34 namespace Utilities
35 {
36  namespace MPI
37  {
63  {
64  public:
68  template <typename NumberType>
69  friend class ::ScaLAPACKMatrix;
70 
78  const unsigned int n_rows,
79  const unsigned int n_columns);
80 
99  const unsigned int n_rows_matrix,
100  const unsigned int n_columns_matrix,
101  const unsigned int row_block_size,
102  const unsigned int column_block_size);
103 
107  ~ProcessGrid();
108 
112  unsigned int
113  get_process_grid_rows() const;
114 
118  unsigned int
119  get_process_grid_columns() const;
120 
126  int
127  get_this_process_row() const;
128 
134  int
135  get_this_process_column() const;
136 
142  template <typename NumberType>
143  void
144  send_to_inactive(NumberType *value, const int count = 1) const;
145 
149  bool
150  is_process_active() const;
151 
152  private:
157  ProcessGrid(MPI_Comm mpi_communicator,
158  const std::pair<unsigned int, unsigned int> &grid_dimensions);
159 
164 
170 
176 
180  const unsigned int this_mpi_process;
181 
185  const unsigned int n_mpi_processes;
186 
191 
196 
203 
210 
215  };
216 
217  /*--------------------- Inline functions --------------------------------*/
218 
219 # ifndef DOXYGEN
220 
221  inline unsigned int
223  {
224  return n_process_rows;
225  }
226 
227 
228 
229  inline unsigned int
231  {
232  return n_process_columns;
233  }
234 
235 
236 
237  inline int
239  {
240  return this_process_row;
241  }
242 
243 
244 
245  inline int
247  {
248  return this_process_column;
249  }
250 
251 
252 
253  inline bool
255  {
256  return mpi_process_is_active;
257  }
258 
259 
260 # endif // ifndef DOXYGEN
261 
262  } // end of namespace MPI
263 
264 } // end of namespace Utilities
265 
266 
267 DEAL_II_NAMESPACE_CLOSE
268 
269 #endif // DEAL_II_WITH_SCALAPACK
270 
271 #endif
unsigned int get_process_grid_rows() const
const unsigned int this_mpi_process
Definition: process_grid.h:180
MPI_Comm mpi_communicator_inactive_with_root
Definition: process_grid.h:169
void send_to_inactive(NumberType *value, const int count=1) const
unsigned int get_process_grid_columns() const
Definition: cuda.h:31
const unsigned int n_mpi_processes
Definition: process_grid.h:185
ProcessGrid(MPI_Comm mpi_communicator, const unsigned int n_rows, const unsigned int n_columns)
int get_this_process_column() const