oomph::LinearAlgebraDistribution Class Reference

#include <linear_algebra_distribution.h>

Public Member Functions

 LinearAlgebraDistribution ()
 
 LinearAlgebraDistribution (const OomphCommunicator &comm, const unsigned &first_row_, const unsigned &n_row_local, const unsigned &n_row=0)
 
 LinearAlgebraDistribution (const OomphCommunicator &comm, const unsigned &n_row, const bool &distributed_=true)
 
 LinearAlgebraDistribution (const OomphCommunicator *const comm_pt, const unsigned &first_row_, const unsigned &n_row_local, const unsigned &n_row=0)
 
 LinearAlgebraDistribution (const OomphCommunicator *const comm_pt, const unsigned &n_row, const bool &distributed_=true)
 
 LinearAlgebraDistribution (const LinearAlgebraDistribution &old_dist)
 Copy Constructor. More...
 
 LinearAlgebraDistribution (const LinearAlgebraDistribution *old_dist_pt)
 pointer based copy constructor More...
 
 ~LinearAlgebraDistribution ()
 Destructor. More...
 
void operator= (const LinearAlgebraDistribution &old_dist)
 Assignment Operator. More...
 
void build (const OomphCommunicator *const comm_pt, const unsigned &first_row, const unsigned &nrow_local, const unsigned &nrow=0)
 
void build (const OomphCommunicator *const comm_pt, const unsigned &nrow, const bool &distributed=true)
 
void build (const LinearAlgebraDistribution &new_dist)
 helper method for the =assignment operator and copy constructor More...
 
void build (const LinearAlgebraDistribution *new_dist_pt)
 
void clear ()
 clears the distribution More...
 
unsigned nrow () const
 access function to the number of global rows. More...
 
unsigned nrow_local () const
 
unsigned nrow_local (const unsigned &p) const
 
unsigned first_row () const
 
unsigned first_row (const unsigned &p) const
 access function for the first row on the p-th processor More...
 
bool distributed () const
 
OomphCommunicatorcommunicator_pt () const
 const access to the communicator pointer More...
 
bool built () const
 
bool operator== (const LinearAlgebraDistribution &other_dist) const
 == Operator More...
 
bool operator!= (const LinearAlgebraDistribution &other_dist) const
 != operator More...
 
unsigned global_to_local_row_map (const unsigned &global_i) const
 return the local index corresponding to the global index More...
 
unsigned rank_of_global_row (const unsigned i) const
 return the processor rank of the global row number i More...
 
Vector< unsignednrow_local_vector () const
 return the nrow_local Vector More...
 
Vector< unsignedfirst_row_vector () const
 return the first_row Vector More...
 

Private Attributes

unsigned Nrow
 the number of global rows More...
 
Vector< unsignedNrow_local
 the number of local rows on the processor More...
 
Vector< unsignedFirst_row
 the first row on this processor More...
 
bool Distributed
 
OomphCommunicatorComm_pt
 the pointer to the MPI communicator object in this distribution More...
 

Friends

std::ostream & operator<< (std::ostream &stream, LinearAlgebraDistribution &dist)
 << operator More...
 

Detailed Description

Describes the distribution of a distributable linear algebra type object. Typically this is a container (such as a DoubleVector) or an operator (e.g Preconditioner or LinearSolver). This object is used in both serial and parallel implementations. In the serial context (no MPI) this just contains an integer indicating the number of rows. In parallel either each processor holds a subset of the set of global rows. (each processor contains only a single continuous block of rows - parametised with variables denoting the first row and the number of local rows) or, all rows are be duplicated across all processors. In parallel this object also contains an OomphCommunicator object which primarily contains the MPI_Comm communicator associated with this object.

Constructor & Destructor Documentation

◆ LinearAlgebraDistribution() [1/7]

oomph::LinearAlgebraDistribution::LinearAlgebraDistribution ( )
inline

Default Constructor - creates a Distribution that has not been setup

68 : Comm_pt(0) {}
OomphCommunicator * Comm_pt
the pointer to the MPI communicator object in this distribution
Definition: linear_algebra_distribution.h:425

◆ LinearAlgebraDistribution() [2/7]

oomph::LinearAlgebraDistribution::LinearAlgebraDistribution ( const OomphCommunicator comm,
const unsigned first_row_,
const unsigned n_row_local,
const unsigned n_row = 0 
)
inline

Constructor. Takes the first_row, nrow_local (both for this processor) and nrow as arguments. If nrow is not provided or equal to 0 then it will be computed automatically

77  : Comm_pt(0)
78  {
79  this->build(&comm, first_row_, n_row_local, n_row);
80  };
void build(const OomphCommunicator *const comm_pt, const unsigned &first_row, const unsigned &nrow_local, const unsigned &nrow=0)
Definition: linear_algebra_distribution.cc:35

References build().

◆ LinearAlgebraDistribution() [3/7]

oomph::LinearAlgebraDistribution::LinearAlgebraDistribution ( const OomphCommunicator comm,
const unsigned n_row,
const bool distributed_ = true 
)
inline

Constructor. Takes the number of global rows and uniformly distributes them over the processors if distributed = true (default), if distributed = false then every row is duplicated on every processor

88  : Comm_pt(0)
89  {
90  this->build(&comm, n_row, distributed_);
91  };

References build().

◆ LinearAlgebraDistribution() [4/7]

oomph::LinearAlgebraDistribution::LinearAlgebraDistribution ( const OomphCommunicator *const  comm_pt,
const unsigned first_row_,
const unsigned n_row_local,
const unsigned n_row = 0 
)
inline

Constructor. Takes the first_row, nrow_local (both for this processor) and nrow as arguments. If nrow is not provided or equal to 0 then it will be computed automatically

100  : Comm_pt(0)
101  {
102  this->build(comm_pt, first_row_, n_row_local, n_row);
103  };

References build().

◆ LinearAlgebraDistribution() [5/7]

oomph::LinearAlgebraDistribution::LinearAlgebraDistribution ( const OomphCommunicator *const  comm_pt,
const unsigned n_row,
const bool distributed_ = true 
)
inline

Constructor. Takes the number of global rows and uniformly distributes them over the processors if distributed = true (default), if distributed = false then every row is duplicated on every processor

111  : Comm_pt(0)
112  {
113  this->build(comm_pt, n_row, distributed_);
114  };

References build().

◆ LinearAlgebraDistribution() [6/7]

oomph::LinearAlgebraDistribution::LinearAlgebraDistribution ( const LinearAlgebraDistribution old_dist)
inline

Copy Constructor.

118  : Comm_pt(0)
119  {
120  this->build(old_dist);
121  }

References build().

◆ LinearAlgebraDistribution() [7/7]

oomph::LinearAlgebraDistribution::LinearAlgebraDistribution ( const LinearAlgebraDistribution old_dist_pt)
inline

pointer based copy constructor

125  : Comm_pt(0)
126  {
127  this->build(old_dist_pt);
128  }

References build().

◆ ~LinearAlgebraDistribution()

oomph::LinearAlgebraDistribution::~LinearAlgebraDistribution ( )
inline

Destructor.

132  {
133  delete Comm_pt;
134  }

References Comm_pt.

Member Function Documentation

◆ build() [1/4]

void oomph::LinearAlgebraDistribution::build ( const LinearAlgebraDistribution new_dist)

helper method for the =assignment operator and copy constructor

Copy the argument distribution. Also a helper method for the =assignment operator and copy constructor

228  {
229  // delete the existing storage
230  First_row.clear();
231  Nrow_local.clear();
232 
233  // if new_dist is not setup
234  if (new_dist.communicator_pt() == 0)
235  {
236  delete Comm_pt;
237  Comm_pt = 0;
238  Distributed = true;
239  Nrow = 0;
240  }
241  else
242  {
243  // copy the communicator
244  delete Comm_pt;
245  Comm_pt = new OomphCommunicator(*new_dist.communicator_pt());
246 
247  // the new distribution is distributed
248  if (new_dist.distributed())
249  {
250  First_row = new_dist.first_row_vector();
251  Nrow_local = new_dist.nrow_local_vector();
252 
253  Distributed = true;
254  }
255  // else if the new ditribution is not distributed
256  else
257  {
258  Distributed = false;
259  }
260  Nrow = new_dist.nrow();
261  }
262  }
Vector< unsigned > Nrow_local
the number of local rows on the processor
Definition: linear_algebra_distribution.h:414
unsigned Nrow
the number of global rows
Definition: linear_algebra_distribution.h:411
Vector< unsigned > First_row
the first row on this processor
Definition: linear_algebra_distribution.h:417
bool Distributed
Definition: linear_algebra_distribution.h:422

References Comm_pt, communicator_pt(), distributed(), Distributed, First_row, first_row_vector(), nrow(), Nrow, Nrow_local, and nrow_local_vector().

◆ build() [2/4]

void oomph::LinearAlgebraDistribution::build ( const LinearAlgebraDistribution new_dist_pt)
inline

Copy the argument distribution. Also a helper method for the =assignment operator and copy constructor

166  {
167  this->build(*new_dist_pt);
168  }

References build().

◆ build() [3/4]

void oomph::LinearAlgebraDistribution::build ( const OomphCommunicator *const  comm_pt,
const unsigned first_row,
const unsigned local_nrow,
const unsigned global_nrow = 0 
)

Sets the distribution. Takes first_row, nrow_local and nrow as arguments. If nrow is not provided or equal to 0 then it is computed automatically

Sets the distribution. Takes first_row, local_nrow and global_nrow as arguments. If global_nrow is not provided or equal to 0 then it is computed automatically

39  {
40  // copy the communicator
41  delete Comm_pt;
42  Comm_pt = new OomphCommunicator(*comm_pt);
43 
44  // get the rank and the number of processors
45  int my_rank = Comm_pt->my_rank();
46  int nproc = Comm_pt->nproc();
47 
48  // resize the storage
49  First_row.clear();
50  First_row.resize(nproc);
51  Nrow_local.clear();
52  Nrow_local.resize(nproc);
53 
54  // set first row and local nrow for this processor
55  First_row[my_rank] = first_row;
56  Nrow_local[my_rank] = local_nrow;
57 
58 #ifdef OOMPH_HAS_MPI
59  // gather the First_row vector
60  unsigned my_nrow_local = Nrow_local[my_rank];
61  MPI_Allgather(&my_nrow_local,
62  1,
63  MPI_UNSIGNED,
64  &Nrow_local[0],
65  1,
66  MPI_UNSIGNED,
67  Comm_pt->mpi_comm());
68 
69  // gather the Nrow_local vector
70  unsigned my_first_row = First_row[my_rank];
71  MPI_Allgather(&my_first_row,
72  1,
73  MPI_UNSIGNED,
74  &First_row[0],
75  1,
76  MPI_UNSIGNED,
77  Comm_pt->mpi_comm());
78 #endif
79 
80  // if global nrow is not provided then compute by summing local_nrow over
81  // all processors
82  if (global_nrow == 0)
83  {
84  if (nproc == 1)
85  {
86  Nrow = local_nrow;
87  }
88  else
89  {
90  Nrow = 0;
91  for (int p = 0; p < nproc; p++)
92  {
93  Nrow += Nrow_local[p];
94  }
95  }
96  }
97  else
98  {
99  Nrow = global_nrow;
100  }
101 
102  // the distribution is true, unless we only have 1 processor
103  if (nproc != 1)
104  {
105  Distributed = true;
106  }
107  else
108  {
109  Distributed = false;
110  }
111 
112 #ifdef OOMPH_HAS_MPI
113 #ifdef PARANOID
114  // paranoid check that the distribution works
115 
116 
117  // check that none of the processors partition overlap
118  for (int p = 0; p < nproc; p++)
119  {
120  if (Nrow_local[p] > 0)
121  {
122  for (int pp = p + 1; pp < nproc; pp++)
123  {
124  if (Nrow_local[pp] > 0)
125  {
126  if ((First_row[p] >= First_row[pp] &&
127  First_row[p] < First_row[pp] + Nrow_local[pp]) ||
128  (First_row[p] + Nrow_local[p] - 1 >= First_row[pp] &&
129  First_row[p] + Nrow_local[p] - 1 <
130  First_row[pp] + Nrow_local[pp]))
131  {
132  std::ostringstream error_message;
133  error_message
134  << "The distributed rows on processor " << p
135  << " and processor " << pp << " overlap.\n"
136  << "Processor " << p << " : first_row = " << First_row[p]
137  << ", nrow = " << Nrow_local[p] << ".\n"
138  << "Processor " << pp << " : first_row = " << First_row[pp]
139  << ", nrow = " << Nrow_local[pp] << ".\n";
140  throw OomphLibWarning(
141  error_message.str(),
142  "LinearAlgebraDistribution::distribute(...)",
144  }
145  }
146  }
147  }
148  }
149 
150  // check that no processor has a row with a global row index greater than
151  // the number of global rows
152  for (int p = 0; p < nproc; p++)
153  {
154  if (First_row[p] + Nrow_local[p] > Nrow)
155  {
156  std::ostringstream error_message;
157  error_message << "Processor " << p << " contains rows " << First_row[p]
158  << " to " << First_row[p] + Nrow_local[p] - 1
159  << " but there are only " << Nrow << " to be distributed."
160  << std::endl;
161  throw OomphLibWarning(error_message.str(),
162  "LinearAlgebraDistribution::distribute(...)",
164  }
165  }
166 #endif
167 #endif
168  }
float * p
Definition: Tutorial_Map_using.cpp:9
unsigned first_row() const
Definition: linear_algebra_distribution.h:261
int my_rank() const
my rank
Definition: communicator.h:176
int nproc() const
number of processors
Definition: communicator.h:157
#define OOMPH_EXCEPTION_LOCATION
Definition: oomph_definitions.h:61

References Comm_pt, Distributed, first_row(), First_row, oomph::OomphCommunicator::my_rank(), oomph::OomphCommunicator::nproc(), Nrow, Nrow_local, OOMPH_EXCEPTION_LOCATION, and p.

Referenced by oomph::PeriodicOrbitAssemblyHandler< NNODE_1D >::adapt_temporal_mesh(), oomph::Problem::assign_eqn_numbers(), build(), oomph::DistributableLinearAlgebraObject::build_distribution(), oomph::LinearAlgebraDistributionHelpers::concatenate(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::FoldHandler::FoldHandler(), oomph::CRDoubleMatrix::get_matrix_transpose(), oomph::HopfHandler::HopfHandler(), LinearAlgebraDistribution(), main(), oomph::SumOfMatrices::multiply(), operator=(), oomph::BlockPitchForkLinearSolver::resolve(), oomph::BlockHopfLinearSolver::solve(), oomph::MumpsSolver::solve(), oomph::FoldHandler::solve_augmented_block_system(), oomph::FoldHandler::solve_block_system(), oomph::HopfHandler::solve_complex_system(), oomph::BlockHopfLinearSolver::solve_for_two_rhs(), oomph::FoldHandler::solve_full_system(), oomph::HopfHandler::solve_full_system(), oomph::HopfHandler::solve_standard_system(), oomph::FoldHandler::~FoldHandler(), and oomph::HopfHandler::~HopfHandler().

◆ build() [4/4]

void oomph::LinearAlgebraDistribution::build ( const OomphCommunicator *const  comm_pt,
const unsigned global_nrow,
const bool distribute = true 
)

Build the LinearAlgebraDistribution. if distributed = true (default) then uniformly distribute nrow over all processors where processors 0 holds approximately the first nrow/n_proc, processor 1 holds the next nrow/n_proc and so on... or if distributed = false then every row is held on every processor

Uniformly distribute global_nrow over all processors where processors 0 holds approximately the first global_nrow/n_proc, processor 1 holds the next global_nrow/n_proc and so on...

179  {
180  // copy the communicator
181  delete Comm_pt;
182  Comm_pt = new OomphCommunicator(*comm_pt);
183 
184  // delete existing storage
185  First_row.clear();
186  Nrow_local.clear();
187 
188  // set global nrow
189  Nrow = global_nrow;
190 
191  // store the distributed flag
192  Distributed = distribute;
193 
194 #ifdef OOMPH_HAS_MPI
195 
196  // if distributed object then compute uniform distribution
197  if (distribute == true)
198  {
199  // get the number of processors
200  int nproc = Comm_pt->nproc();
201 
202  // resize the storage
203  First_row.resize(nproc);
204  Nrow_local.resize(nproc);
205 
206  // compute first row
207  for (int p = 0; p < nproc; p++)
208  {
209  First_row[p] =
210  unsigned((double(p) / double(nproc)) * double(global_nrow));
211  }
212 
213  // compute local nrow
214  for (int p = 0; p < nproc - 1; p++)
215  {
216  Nrow_local[p] = First_row[p + 1] - First_row[p];
217  }
218  Nrow_local[nproc - 1] = global_nrow - First_row[nproc - 1];
219  }
220 #endif
221  }

References Comm_pt, Distributed, First_row, oomph::OomphCommunicator::nproc(), Nrow, Nrow_local, and p.

◆ built()

◆ clear()

void oomph::LinearAlgebraDistribution::clear ( )
inline

clears the distribution

172  {
173  // delete the communicator
174  delete Comm_pt;
175  Comm_pt = 0;
176 
177  // delete first_row and nrow_local
178  First_row.clear();
179  Nrow_local.clear();
180 
181  // zero Nrow
182  Nrow = 0;
183  }

References Comm_pt, First_row, Nrow, and Nrow_local.

Referenced by oomph::DistributableLinearAlgebraObject::clear_distribution(), oomph::DoubleVectorHelpers::concatenate_without_communication(), oomph::MGSolver< DIM >::setup_mg_structures(), and oomph::DoubleVectorHelpers::split_without_communication().

◆ communicator_pt()

OomphCommunicator* oomph::LinearAlgebraDistribution::communicator_pt ( ) const
inline

const access to the communicator pointer

336  {
337  return Comm_pt;
338  }

References Comm_pt.

Referenced by build(), oomph::Problem::calculate_continuation_derivatives_fd_helper(), oomph::Problem::calculate_continuation_derivatives_helper(), oomph::CRDoubleMatrixHelpers::concatenate(), oomph::DoubleVectorHelpers::concatenate(), oomph::DoubleVectorHelpers::concatenate_without_communication(), oomph::CRDoubleMatrixHelpers::concatenate_without_communication(), oomph::TrilinosEpetraHelpers::copy_to_oomphlib_vector(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_vector(), oomph::TrilinosEpetraHelpers::create_epetra_map(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::HypreHelpers::create_HYPRE_Vector(), oomph::DoubleVectorHaloScheme::DoubleVectorHaloScheme(), oomph::SuperLUSolver::factorise(), oomph::MumpsSolver::factorise(), oomph::Problem::get_jacobian(), oomph::CRDoubleMatrix::global_matrix(), oomph::HypreInterface::hypre_matrix_setup(), oomph::HypreInterface::hypre_solver_setup(), oomph::BlockPreconditioner< MATRIX >::internal_get_block(), oomph::MGSolver< DIM >::interpolation_matrix_set(), oomph::HelmholtzMGPreconditioner< DIM >::interpolation_matrix_set(), oomph::TrilinosEpetraHelpers::multiply(), oomph::DenseDoubleMatrix::multiply(), oomph::CCDoubleMatrix::multiply(), oomph::DenseDoubleMatrix::multiply_transpose(), oomph::CCDoubleMatrix::multiply_transpose(), oomph::Problem::newton_solve_continuation(), operator==(), oomph::DoubleMultiVector::output(), oomph::DoubleVector::output(), oomph::DoubleMultiVector::redistribute(), oomph::DoubleVector::redistribute(), oomph::CRDoubleMatrix::redistribute(), oomph::MGSolver< DIM >::self_test(), oomph::MatrixVectorProduct::setup(), oomph::Preconditioner::setup(), oomph::MGSolver< DIM >::setup_mg_structures(), oomph::HelmholtzMGPreconditioner< DIM >::setup_mg_structures(), oomph::MGSolver< DIM >::setup_smoothers(), oomph::HelmholtzMGPreconditioner< DIM >::setup_smoothers(), oomph::DenseLU::solve(), oomph::TrilinosAztecOOSolver::solve(), oomph::MumpsSolver::solve(), oomph::AugmentedProblemGMRES::solve_helper(), oomph::DoubleVectorHelpers::split(), and oomph::DoubleVectorHelpers::split_without_communication().

◆ distributed()

bool oomph::LinearAlgebraDistribution::distributed ( ) const
inline

access function to the distributed - indicates whether the distribution is serial or distributed

330  {
331  return Distributed;
332  }

References Distributed.

Referenced by oomph::SuperLUSolver::backsub_serial(), oomph::SuperLUSolver::backsub_transpose_serial(), build(), oomph::Problem::calculate_continuation_derivatives_fd_helper(), oomph::Problem::calculate_continuation_derivatives_helper(), oomph::Smoother::check_validity_of_solve_helper_inputs(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_vector(), oomph::TrilinosEpetraHelpers::create_epetra_map(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::HypreHelpers::create_HYPRE_Vector(), oomph::DistributableLinearAlgebraObject::distributed(), oomph::DoubleVectorHaloScheme::DoubleVectorHaloScheme(), oomph::Problem::get_eigenproblem_matrices(), oomph::Problem::get_jacobian(), oomph::BlockPreconditioner< MATRIX >::internal_get_block(), main(), oomph::Problem::newton_solve_continuation(), operator==(), oomph::DoubleMultiVector::redistribute(), oomph::DoubleVector::redistribute(), oomph::CRDoubleMatrix::redistribute(), oomph::MatrixVectorProduct::setup(), oomph::DenseLU::solve(), oomph::SuperLUSolver::solve(), oomph::MumpsSolver::solve(), oomph::CG< MATRIX >::solve_helper(), oomph::BiCGStab< MATRIX >::solve_helper(), oomph::GMRES< MATRIX >::solve_helper(), oomph::AugmentedProblemGMRES::solve_helper(), oomph::HelmholtzGMRESMG< MATRIX >::solve_helper(), oomph::HelmholtzFGMRESMG< MATRIX >::solve_helper(), and oomph::SuperLUSolver::solve_transpose().

◆ first_row() [1/2]

unsigned oomph::LinearAlgebraDistribution::first_row ( ) const
inline

access function for the first row on this processor. If not distributed then this is just zero.

262  {
263  // return the first row
264 #ifdef OOMPH_HAS_MPI
265  if (Distributed)
266  {
267 #ifdef PARANOID
268  if (Comm_pt == 0)
269  {
270  throw OomphLibError(
271  "LinearAlgebraDistribution has not been built : Comm_pt == 0.",
274  }
275 #endif
276 
277  return First_row[Comm_pt->my_rank()];
278  }
279  else
280  {
281  return 0;
282  }
283 #else
284  return 0;
285 #endif
286  }
#define OOMPH_CURRENT_FUNCTION
Definition: oomph_definitions.h:86

References Comm_pt, Distributed, First_row, oomph::OomphCommunicator::my_rank(), OOMPH_CURRENT_FUNCTION, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::NavierStokesSchurComplementPreconditioner::assemble_inv_press_and_veloc_mass_matrix_diagonal(), oomph::PressureBasedSolidLSCPreconditioner::assemble_mass_matrix_diagonal(), build(), oomph::CRDoubleMatrixHelpers::concatenate(), oomph::DoubleVectorHelpers::concatenate(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_vector(), oomph::TrilinosEpetraHelpers::create_epetra_map(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::HypreHelpers::create_HYPRE_Vector(), create_matrix_ascend_col_row(), create_vector_ascend_row(), oomph::DoubleVectorHaloScheme::DoubleVectorHaloScheme(), oomph::DistributableLinearAlgebraObject::first_row(), oomph::Problem::global_dof_pt(), oomph::HypreInterface::hypre_solve(), oomph::BlockPreconditioner< MATRIX >::index_in_block(), oomph::BlockPreconditioner< MATRIX >::internal_dof_number(), oomph::BlockPreconditioner< MATRIX >::internal_index_in_dof(), main(), operator==(), oomph::PitchForkHandler::PitchForkHandler(), rank_of_global_row(), oomph::DoubleMultiVector::redistribute(), oomph::DoubleVector::redistribute(), oomph::CRDoubleMatrix::redistribute(), and oomph::CRDoubleMatrix::sparse_indexed_output_with_offset().

◆ first_row() [2/2]

unsigned oomph::LinearAlgebraDistribution::first_row ( const unsigned p) const
inline

access function for the first row on the p-th processor

290  {
291 #ifdef PARANOID
292  if (Comm_pt == 0)
293  {
294  throw OomphLibError(
295  "LinearAlgebraDistribution has not been built : Comm_pt == 0.",
298  }
299  if (p >= unsigned(Comm_pt->nproc()))
300  {
301  std::ostringstream error_message;
302  error_message << "Requested first_row(" << p
303  << "), but this distribution is defined "
304  << "on " << Comm_pt->nproc() << "processors.";
305  throw OomphLibError(error_message.str(),
308  }
309 
310 #endif
311 
312  // return the first row
313 #ifdef OOMPH_HAS_MPI
314  if (Distributed)
315  {
316  return First_row[p];
317  }
318  else
319  {
320  return 0;
321  }
322 #else
323  return 0;
324 #endif
325  }

References Comm_pt, Distributed, First_row, oomph::OomphCommunicator::nproc(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and p.

◆ first_row_vector()

Vector<unsigned> oomph::LinearAlgebraDistribution::first_row_vector ( ) const
inline

return the first_row Vector

405  {
406  return First_row;
407  }

References First_row.

Referenced by build().

◆ global_to_local_row_map()

unsigned oomph::LinearAlgebraDistribution::global_to_local_row_map ( const unsigned global_i) const
inline

return the local index corresponding to the global index

366  {
367 #ifdef PARANOID
368  if (global_i >= Nrow)
369  {
370  throw OomphLibError(
371  "Requested global row outside the number of global rows",
374  }
375 #endif
376  int local_i = static_cast<int>(global_i);
377  int p = 0;
378  while ((int)(local_i - (int)nrow_local(p)) >= 0)
379  {
380  local_i -= (int)nrow_local(p);
381  p++;
382  }
383  return (unsigned)local_i;
384  }
unsigned nrow_local() const
Definition: linear_algebra_distribution.h:193
return int(ret)+1

References int(), Nrow, nrow_local(), OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and p.

◆ nrow()

◆ nrow_local() [1/2]

unsigned oomph::LinearAlgebraDistribution::nrow_local ( ) const
inline

access function for the num of local rows on this processor. If no MPI then Nrow is returned.

194  {
195  // return the nrow_local
196 #ifdef OOMPH_HAS_MPI
197  if (Distributed)
198  {
199 #ifdef PARANOID
200  if (Comm_pt == 0)
201  {
202  throw OomphLibError(
203  "LinearAlgebraDistribution has not been built : Comm_pt == 0.",
204  "LinearAlgebraDistribution::nrow_local()",
206  }
207 #endif
208 
209  return Nrow_local[Comm_pt->my_rank()];
210  }
211  else
212  {
213  return Nrow;
214  }
215 #else
216  return Nrow;
217 #endif
218  }

References Comm_pt, Distributed, oomph::OomphCommunicator::my_rank(), Nrow, Nrow_local, and OOMPH_EXCEPTION_LOCATION.

Referenced by oomph::Problem::adapt(), oomph::Problem::adaptive_unsteady_newton_solve(), oomph::OomphLibPreconditionerEpetraOperator::ApplyInverse(), oomph::Problem::arc_length_step_solve_helper(), oomph::NavierStokesSchurComplementPreconditioner::assemble_inv_press_and_veloc_mass_matrix_diagonal(), oomph::PressureBasedSolidLSCPreconditioner::assemble_mass_matrix_diagonal(), oomph::Problem::calculate_continuation_derivatives_fd_helper(), oomph::Problem::calculate_continuation_derivatives_helper(), oomph::CRDoubleMatrixHelpers::concatenate(), oomph::CRDoubleMatrixHelpers::concatenate_without_communication(), oomph::CRDoubleMatrix::CRDoubleMatrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix(), oomph::TrilinosEpetraHelpers::create_distributed_epetra_matrix_for_aztecoo(), oomph::TrilinosEpetraHelpers::create_epetra_map(), oomph::HypreHelpers::create_HYPRE_Matrix(), oomph::HypreHelpers::create_HYPRE_Vector(), create_matrix_ascend_col_row(), create_vector_ascend_row(), oomph::Problem::get_eigenproblem_matrices(), oomph::Problem::get_hessian_vector_products(), oomph::Problem::get_jacobian(), oomph::Problem::get_residuals(), oomph::Problem::global_dof_pt(), global_to_local_row_map(), oomph::HypreInterface::hypre_solve(), oomph::BlockPreconditioner< MATRIX >::index_in_block(), oomph::BlockPreconditioner< MATRIX >::internal_dof_number(), oomph::BlockPreconditioner< MATRIX >::internal_index_in_dof(), main(), oomph::Problem::newton_solve(), oomph::Problem::newton_solve_continuation(), oomph::DistributableLinearAlgebraObject::nrow_local(), operator==(), oomph::PitchForkHandler::PitchForkHandler(), rank_of_global_row(), oomph::DoubleMultiVector::redistribute(), oomph::DoubleVector::redistribute(), oomph::CRDoubleMatrix::redistribute(), oomph::BlockPitchForkLinearSolver::resolve(), oomph::Problem::restore_dof_values(), oomph::TwoPlusThreeUpperTriangularWithReplace< MATRIX >::setup(), oomph::CoarseTwoPlusTwoPlusOne< MATRIX >::setup(), oomph::OnePlusFourWithTwoCoarse< MATRIX >::setup(), oomph::BlockPitchForkLinearSolver::solve(), oomph::PitchForkHandler::solve_block_system(), oomph::PitchForkHandler::solve_full_system(), oomph::Problem::store_current_dof_values(), and oomph::PitchForkHandler::~PitchForkHandler().

◆ nrow_local() [2/2]

unsigned oomph::LinearAlgebraDistribution::nrow_local ( const unsigned p) const
inline

access function for the num of local rows on this processor. If no MPI the nrow is returned

223  {
224 #ifdef PARANOID
225  if (Comm_pt == 0)
226  {
227  throw OomphLibError(
228  "LinearAlgebraDistribution has not been built : Comm_pt == 0.",
231  }
232  if (p >= unsigned(Comm_pt->nproc()))
233  {
234  std::ostringstream error_message;
235  error_message << "Requested nrow_local(" << p
236  << "), but this distribution is defined "
237  << "on " << Comm_pt->nproc() << "processors.";
238  throw OomphLibError(error_message.str(),
241  }
242 #endif
243 
244  // return the nrow_local
245 #ifdef OOMPH_HAS_MPI
246  if (Distributed)
247  {
248  return Nrow_local[p];
249  }
250  else
251  {
252  return Nrow;
253  }
254 #else
255  return Nrow;
256 #endif
257  }

References Comm_pt, Distributed, oomph::OomphCommunicator::nproc(), Nrow, Nrow_local, OOMPH_CURRENT_FUNCTION, OOMPH_EXCEPTION_LOCATION, and p.

◆ nrow_local_vector()

Vector<unsigned> oomph::LinearAlgebraDistribution::nrow_local_vector ( ) const
inline

return the nrow_local Vector

399  {
400  return Nrow_local;
401  }

References Nrow_local.

Referenced by build().

◆ operator!=()

bool oomph::LinearAlgebraDistribution::operator!= ( const LinearAlgebraDistribution other_dist) const
inline

!= operator

356  {
357  return !(*this == other_dist);
358  }

◆ operator=()

void oomph::LinearAlgebraDistribution::operator= ( const LinearAlgebraDistribution old_dist)
inline

Assignment Operator.

138  {
139  this->build(old_dist);
140  }

References build().

◆ operator==()

bool oomph::LinearAlgebraDistribution::operator== ( const LinearAlgebraDistribution other_dist) const

== Operator

operator==

270  {
271 #ifdef OOMPH_HAS_MPI
272  // compare the communicators
273  if (!((*Comm_pt) == (*other_dist.communicator_pt())))
274  {
275  return false;
276  }
277 
278  // compare Distributed
279  if (Distributed != other_dist.distributed())
280  {
281  return false;
282  }
283 
284  // if not distributed compare nrow
285  if (!Distributed)
286  {
287  if (other_dist.nrow() == Nrow)
288  {
289  return true;
290  }
291  return false;
292  }
293 
294  // compare
295  bool flag = true;
296  int nproc = Comm_pt->nproc();
297  for (int i = 0; i < nproc && flag == true; i++)
298  {
299  if (other_dist.first_row(i) != First_row[i] ||
300  other_dist.nrow_local(i) != Nrow_local[i])
301  {
302  flag = false;
303  }
304  }
305  return flag;
306 #else
307  if (other_dist.nrow() == Nrow)
308  {
309  return true;
310  }
311  return false;
312 #endif
313  }
int i
Definition: BiCGSTAB_step_by_step.cpp:9

References Comm_pt, communicator_pt(), distributed(), Distributed, first_row(), First_row, i, oomph::OomphCommunicator::nproc(), nrow(), Nrow, nrow_local(), and Nrow_local.

◆ rank_of_global_row()

unsigned oomph::LinearAlgebraDistribution::rank_of_global_row ( const unsigned  i) const
inline

return the processor rank of the global row number i

388  {
389  unsigned p = 0;
390  while (i < first_row(p) || i >= first_row(p) + nrow_local(p))
391  {
392  p++;
393  }
394  return p;
395  }

References first_row(), i, nrow_local(), and p.

Referenced by oomph::CRDoubleMatrixHelpers::concatenate(), oomph::DoubleVectorHelpers::concatenate(), oomph::DoubleVectorHaloScheme::DoubleVectorHaloScheme(), and oomph::BlockPreconditioner< MATRIX >::index_in_block().

Friends And Related Function Documentation

◆ operator<<

std::ostream& operator<< ( std::ostream &  stream,
LinearAlgebraDistribution dist 
)
friend

<< operator

320  {
321  stream << "nrow()=" << dist.nrow() << ", first_row()=" << dist.first_row()
322  << ", nrow_local()=" << dist.nrow_local()
323  << ", distributed()=" << dist.distributed() << std::endl;
324  return stream;
325  }

Member Data Documentation

◆ Comm_pt

OomphCommunicator* oomph::LinearAlgebraDistribution::Comm_pt
private

the pointer to the MPI communicator object in this distribution

Referenced by build(), built(), clear(), communicator_pt(), first_row(), nrow_local(), operator==(), and ~LinearAlgebraDistribution().

◆ Distributed

bool oomph::LinearAlgebraDistribution::Distributed
private

flag to indicate whether this distribution describes an object that is distributed over the processors of Comm_pt (true) or duplicated over the processors of Comm_pt (false)

Referenced by build(), distributed(), first_row(), nrow_local(), and operator==().

◆ First_row

Vector<unsigned> oomph::LinearAlgebraDistribution::First_row
private

the first row on this processor

Referenced by build(), clear(), first_row(), first_row_vector(), and operator==().

◆ Nrow

unsigned oomph::LinearAlgebraDistribution::Nrow
private

the number of global rows

Referenced by build(), clear(), global_to_local_row_map(), nrow(), nrow_local(), and operator==().

◆ Nrow_local

Vector<unsigned> oomph::LinearAlgebraDistribution::Nrow_local
private

the number of local rows on the processor

Referenced by build(), clear(), nrow_local(), nrow_local_vector(), and operator==().


The documentation for this class was generated from the following files: