communicator.h
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-2025 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#ifndef OOMPH_COMMUNICATOR_CLASS_HEADER
27#define OOMPH_COMMUNICATOR_CLASS_HEADER
28
29// Config header
30#ifdef HAVE_CONFIG_H
31#include <oomph-lib-config.h>
32#endif
33
34// MPI headers
35#ifdef OOMPH_HAS_MPI
36#include "mpi.h"
37#endif
38
39// Oomph-lib error handler
40#include "oomph_definitions.h"
41#include "oomph_utilities.h"
42
43namespace oomph
44{
45 template<class T>
46 class DenseMatrix;
47
48 //=========================================================================
49 /// An oomph-lib wrapper to the MPI_Comm communicator object. Just
50 /// contains an MPI_Comm object (which is a pointer) and wrappers to
51 /// the MPI_... methods.
52 //=========================================================================
54 {
55 public:
56#ifdef OOMPH_HAS_MPI
57 /// Construct a communicator from a MPI_Comm object.
58 /// if the bool owns_mpi_comm is true then this communicator is responsible
59 /// for the destruction of the mpi_communicator. The mpi comm will be freed
60 /// when the destructor is called.
62 const bool& owns_mpi_comm = false)
64 {
65 // store a pointer to the communicator
67
68 // hold owns_mpi_comm
70 }
71#endif
72
73 /// Serial constructor
75#ifdef OOMPH_HAS_MPI
77#endif
78 {
79 }
80
81 /// Copy constructor.
83#ifdef OOMPH_HAS_MPI
84 : Owns_mpi_comm(false)
85 {
86 if (communicator.serial_communicator())
87 {
89 }
90 else
91 {
92 Comm = communicator.mpi_comm();
93 Serial_communicator = false;
94 }
95 }
96#else
97 {
98 }
99#endif
100
101 /// Pointer (copy) constructor.
102 OomphCommunicator(const OomphCommunicator* communicator_pt)
103#ifdef OOMPH_HAS_MPI
104 : Owns_mpi_comm(false)
105 {
106 if (communicator_pt->serial_communicator())
107 {
108 Serial_communicator = true;
109 }
110 else
111 {
112 Comm = communicator_pt->mpi_comm();
113 Serial_communicator = false;
114 }
115 }
116#else
117 {
118 }
119#endif
120
121 /// Destructor. If MPI and this preconditioner owns the MPI_Comm
122 /// object then MPI_Comm_free is called, otherwise nothing happens as
123 /// the destruction of the underlying MPI_Comm object is the responsibility
124 /// of another communicator.
125 ~OomphCommunicator()
126 {
127#ifdef OOMPH_HAS_MPI
128 if (Owns_mpi_comm)
129 {
130 MPI_Comm_free(&Comm);
131 }
132#endif
133 }
134
135 /// assignment operator
136 void operator=(const OomphCommunicator& communicator)
137 {
138#ifdef OOMPH_HAS_MPI
139 if (Owns_mpi_comm)
140 {
141 MPI_Comm_free(&Comm);
142 }
143 Owns_mpi_comm = false;
144 if (communicator.serial_communicator())
145 {
146 Serial_communicator = true;
147 }
148 else
149 {
150 Serial_communicator = false;
151 Comm = communicator.mpi_comm();
152 }
153#endif
154 }
155
156 /// number of processors
157 int nproc() const
158 {
159#ifdef OOMPH_HAS_MPI
161 {
162 return 1;
163 }
164 else
165 {
166 int n_proc = 1;
167 MPI_Comm_size(Comm, &n_proc);
168 return n_proc;
169 }
170#else
171 return 1;
172#endif
173 }
174
175 /// my rank
176 int my_rank() const
177 {
178#ifdef OOMPH_HAS_MPI
180 {
181 return 0;
182 }
183 else
184 {
185 int My_rank = 0;
186 MPI_Comm_rank(Comm, &My_rank);
187 return My_rank;
188 }
189#else
190 return 0;
191#endif
192 }
193
194 /// == operator - only returns true if communicators are MPI_IDENT,
195 /// i.e. if both group and context are the same
196 bool operator==(const OomphCommunicator& other_comm) const
197 {
198#ifdef OOMPH_HAS_MPI
199 if (Serial_communicator != other_comm.serial_communicator())
200 {
201 return false;
202 }
203 else if (Serial_communicator)
204 {
205 return true;
206 }
207 else
208 {
209 int flag;
210 MPI_Comm_compare(Comm, other_comm.mpi_comm(), &flag);
211 if (flag == MPI_IDENT)
212 {
213 return true;
214 }
215 return false;
216 }
217#else
218 return true;
219#endif
220 }
221
222 /// != operator returns !(==operator) (see ==operator for more
223 /// details)
224 bool operator!=(const OomphCommunicator& other_comm) const
225 {
226 return !(*this == other_comm);
227 }
228
229#ifdef OOMPH_HAS_MPI
230 /// split the communicator: color is an integer label for the sub
231 /// group of processors. key is the rank in this sub group
232 OomphCommunicator* split(const int& color, const int& key)
233 {
234#ifdef PARANOID
236 {
237 std::ostringstream error_message_stream;
238 error_message_stream << "Attempted to split a serial communicator.";
239 throw OomphLibError(error_message_stream.str(),
240 OOMPH_CURRENT_FUNCTION,
241 OOMPH_EXCEPTION_LOCATION);
242 }
243#endif
244
245 // pointer for the split oomph-lib communicator
246 OomphCommunicator* split_comm_pt;
247
248 // the pointer for the new mpi communicator
249 MPI_Comm* mpi_comm_pt = new MPI_Comm;
250
251 // get the split communicator
252 MPI_Comm_split(Comm, color, key, mpi_comm_pt);
253
254 // assemble the new oomph-lib ocmmunicator
255 // the new oomph-lib communicator owns the MPI_Comm at mpi_comm_pt
256 // and is responsible for its destruction
257 split_comm_pt = new OomphCommunicator(*mpi_comm_pt, true);
258
259 // and return
260 return split_comm_pt;
261 }
262#endif
263
264
265#ifdef OOMPH_HAS_MPI
266 /// access function to the underlying MPI_Comm object
267 MPI_Comm mpi_comm() const
268 {
269#ifdef PARANOID
271 {
272 std::ostringstream error_message_stream;
273 error_message_stream
274 << "Requested the MPI_Comm object for a serial communicator.";
275 throw OomphLibError(error_message_stream.str(),
276 OOMPH_CURRENT_FUNCTION,
277 OOMPH_EXCEPTION_LOCATION);
278 }
279#endif
280 return Comm;
281 }
282
283
284 /// Access function (read only) to the serial communicator flag
285 bool serial_communicator() const
286 {
287 return Serial_communicator;
288 }
289
290 /// broadcast function for Vector<int>
291 void broadcast(const int& source, Vector<int>& x)
292 {
293 int n;
294
295 // Get number of entries on processor source (where the vector exists)
296 unsigned long n_long;
297 if (this->my_rank() == source)
298 {
299 n_long = x.size();
300
301 // Convert to int
302 n = static_cast<int>(n_long);
303 }
304
305 // Broadcast to everybody how many entries to expect
306 MPI_Bcast(&n, 1, MPI_INT, source, this->mpi_comm());
307
308 // Resize Vector everywhere else in preparation
309 if (this->my_rank() != source)
310 {
311 x.resize(n);
312 }
313
314 // Broadcast the Vector directly
315 MPI_Bcast(&x[0], n, MPI_INT, source, this->mpi_comm());
316 }
317
318 /// broadcast function for Vector<double>
319 void broadcast(const int& source, Vector<double>& x)
320 {
321 int n;
322
323 // Get number of entries on processor source (where the vector exists)
324 unsigned long n_long;
325 if (this->my_rank() == source)
326 {
327 n_long = x.size();
328
329 // Convert to int
330 n = static_cast<int>(n_long);
331 }
332
333 // Broadcast to everybody how many entries to expect
334 MPI_Bcast(&n, 1, MPI_INT, source, this->mpi_comm());
335
336 // Resize Vector everywhere else in preparation
337 if (this->my_rank() != source)
338 {
339 x.resize(n);
340 }
341
342 // Broadcast the Vector directly
343 MPI_Bcast(&x[0], n, MPI_DOUBLE, source, this->mpi_comm());
344 }
345
346 /// broadcast function for DenseMatrix<double>
347 void broadcast(const int& source, DenseMatrix<double>& x);
348
349#endif
350
351 private:
352#ifdef OOMPH_HAS_MPI
353
354 /// the MPI_Comm communicator
356
357 /// boolean indiacting whether this communicator owns the underlying
358 /// MPI_Comm - if so the destructor will free
359 /// the communicator
361
362 /// boolean to indicate if this communicator is for serial problems.
363 /// This is true when serial codes are compiled under MPI
365
366#endif
367
368 }; // end of class Communicator
369
370} // namespace oomph
371
372#endif
double size() const
Calculate the size of the element (length, area, volume,...) in Eulerian computational coordinates....
Definition elements.cc:4320
An oomph-lib wrapper to the MPI_Comm communicator object. Just contains an MPI_Comm object (which is ...
MPI_Comm Comm
the MPI_Comm communicator
OomphCommunicator()
Serial constructor.
bool Owns_mpi_comm
boolean indiacting whether this communicator owns the underlying MPI_Comm - if so the destructor will...
bool Serial_communicator
boolean to indicate if this communicator is for serial problems. This is true when serial codes are c...
OomphCommunicator(const MPI_Comm mpi_communicator, const bool &owns_mpi_comm=false)
Construct a communicator from a MPI_Comm object. if the bool owns_mpi_comm is true then this communic...
TAdvectionDiffusionReactionElement<NREAGENT,DIM,NNODE_1D> elements are isoparametric triangular DIM-d...
void split(const DoubleVector &in_vector, Vector< DoubleVector * > &out_vector_pt)
Split a DoubleVector into the out DoubleVectors. Let vec_A be the in Vector, and let vec_B and vec_C ...
DRAIG: Change all instances of (SPATIAL_DIM) to (DIM-1).