communicator.cc
Go to the documentation of this file.
1 // LIC// ====================================================================
2 // LIC// This file forms part of oomph-lib, the object-oriented,
3 // LIC// multi-physics finite-element library, available
4 // LIC// at http://www.oomph-lib.org.
5 // LIC//
6 // LIC// Copyright (C) 2006-2023 Matthias Heil and Andrew Hazel
7 // LIC//
8 // LIC// This library is free software; you can redistribute it and/or
9 // LIC// modify it under the terms of the GNU Lesser General Public
10 // LIC// License as published by the Free Software Foundation; either
11 // LIC// version 2.1 of the License, or (at your option) any later version.
12 // LIC//
13 // LIC// This library is distributed in the hope that it will be useful,
14 // LIC// but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // LIC// MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // LIC// Lesser General Public License for more details.
17 // LIC//
18 // LIC// You should have received a copy of the GNU Lesser General Public
19 // LIC// License along with this library; if not, write to the Free Software
20 // LIC// Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA
21 // LIC// 02110-1301 USA.
22 // LIC//
23 // LIC// The authors may be contacted at oomph-lib@maths.man.ac.uk.
24 // LIC//
25 // LIC//====================================================================
26 // Functions for OomphCommunicator
27 
28 // MPI headers
29 #ifdef OOMPH_HAS_MPI
30 #include "mpi.h"
31 #endif
32 
33 // Oomph-lib error handler
34 #include "communicator.h"
35 #include "matrices.h"
36 
37 namespace oomph
38 {
39 #ifdef OOMPH_HAS_MPI
40 
41  //=============================================================================
42  /// A broadcast function for DenseMatrix<double>
43  //=============================================================================
44  void OomphCommunicator::broadcast(const int& source, DenseMatrix<double>& x)
45  {
46  // Get number of entries on processor source (where the matrix exists)
47  unsigned nrow, ncol;
48  if (this->my_rank() == source)
49  {
50  nrow = x.nrow();
51  if (nrow > 0)
52  {
53  ncol = x.ncol();
54  }
55  else
56  {
57  ncol = 0;
58  }
59  }
60 
61  // Broadcast to everybody how many entries to expect
62  MPI_Bcast(&nrow, 1, MPI_UNSIGNED_LONG, source, this->mpi_comm());
63  MPI_Bcast(&ncol, 1, MPI_UNSIGNED_LONG, source, this->mpi_comm());
64 
65  if (ncol != 0 && nrow != 0)
66  {
67  // convert into a C-style array
68  double* x_bcast = new double[nrow * ncol];
69 
70  if (this->my_rank() == source)
71  for (unsigned long i = 0; i < nrow; i++)
72  {
73  for (unsigned long j = 0; j < ncol; j++)
74  {
75  x_bcast[i * ncol + j] = x(i, j);
76  }
77  }
78 
79  // broadcast the array
80  MPI_Bcast(x_bcast, ncol * nrow, MPI_DOUBLE, source, this->mpi_comm());
81 
82  // Now convert back into matrix (everywhere apart from source)
83  if (this->my_rank() != source)
84  {
85  x.resize(nrow, ncol);
86  for (unsigned long i = 0; i < nrow; i++)
87  for (unsigned long j = 0; j < ncol; j++)
88  {
89  x(i, j) = x_bcast[i * ncol + j];
90  }
91  }
92  // delete C-style array
93  delete[] x_bcast;
94  }
95  else
96  {
97  x.resize(0, 0);
98  }
99  }
100 
101 #endif
102 
103 } // namespace oomph
cstr elem_len * i
Definition: cfortran.h:603
//////////////////////////////////////////////////////////////////// ////////////////////////////////...