communicator.C
Go to the documentation of this file.
1 // The libMesh Finite Element Library.
2 // Copyright (C) 2002-2018 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner
3 
4 // This library is free software; you can redistribute it and/or
5 // modify it under the terms of the GNU Lesser General Public
6 // License as published by the Free Software Foundation; either
7 // version 2.1 of the License, or (at your option) any later version.
8 
9 // This library is distributed in the hope that it will be useful,
10 // but WITHOUT ANY WARRANTY; without even the implied warranty of
11 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 // Lesser General Public License for more details.
13 
14 // You should have received a copy of the GNU Lesser General Public
15 // License along with this library; if not, write to the Free Software
16 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
17 
18 
19 // Parallel includes
20 #include "libmesh/communicator.h"
21 
22 // libMesh includes
24 
25 namespace libMesh
26 {
27 
28 namespace Parallel
29 {
30 
31 
32 // ------------------------------------------------------------
33 // Simple Communicator member functions
34 
35 void Communicator::reference_unique_tag(int tagvalue) const
36 {
37  // This had better be an already-acquired tag.
38  libmesh_assert(used_tag_values.count(tagvalue));
39 
40  used_tag_values[tagvalue]++;
41 }
42 
43 
44 void Communicator::dereference_unique_tag(int tagvalue) const
45 {
46  // This had better be an already-acquired tag.
47  libmesh_assert(used_tag_values.count(tagvalue));
48 
49  used_tag_values[tagvalue]--;
50  // If we don't have any more outstanding references, we
51  // don't even need to keep this tag in our "used" set.
52  if (!used_tag_values[tagvalue])
53  used_tag_values.erase(tagvalue);
54 }
55 
56 
58 #ifdef LIBMESH_HAVE_MPI
59  _communicator(MPI_COMM_SELF),
60 #endif
61  _rank(0),
62  _size(1),
63  _send_mode(DEFAULT),
64  used_tag_values(),
65  _I_duped_it(false) {}
66 
67 
69 #ifdef LIBMESH_HAVE_MPI
70  _communicator(MPI_COMM_SELF),
71 #endif
72  _rank(0),
73  _size(1),
74  _send_mode(DEFAULT),
75  used_tag_values(),
76  _I_duped_it(false)
77 {
78  this->assign(comm);
79 }
80 
81 
83 {
84  this->clear();
85 }
86 
87 
88 #ifdef LIBMESH_HAVE_MPI
89 void Communicator::split(int color, int key, Communicator & target) const
90 {
91  target.clear();
92  MPI_Comm newcomm;
93  libmesh_call_mpi
94  (MPI_Comm_split(this->get(), color, key, &newcomm));
95 
96  target.assign(newcomm);
97  target._I_duped_it = (color != MPI_UNDEFINED);
98  target.send_mode(this->send_mode());
99 }
100 #else
101 void Communicator::split(int, int, Communicator & target) const
102 {
103  target.assign(this->get());
104 }
105 #endif
106 
107 
109 {
110  this->duplicate(comm._communicator);
111  this->send_mode(comm.send_mode());
112 }
113 
114 
115 #ifdef LIBMESH_HAVE_MPI
117 {
118  if (_communicator != MPI_COMM_NULL)
119  {
120  libmesh_call_mpi
121  (MPI_Comm_dup(comm, &_communicator));
122 
123  _I_duped_it = true;
124  }
125  this->assign(_communicator);
126 }
127 #else
128 void Communicator::duplicate(const communicator &) { }
129 #endif
130 
131 
133 #ifdef LIBMESH_HAVE_MPI
134  if (_I_duped_it)
135  {
136  libmesh_assert (_communicator != MPI_COMM_NULL);
137  libmesh_call_mpi
138  (MPI_Comm_free(&_communicator));
139 
140  _communicator = MPI_COMM_NULL;
141  }
142  _I_duped_it = false;
143 #endif
144 }
145 
146 
148 {
149  this->clear();
150  this->assign(comm);
151  return *this;
152 }
153 
154 
156 {
157  _communicator = comm;
158 #ifdef LIBMESH_HAVE_MPI
159  if (_communicator != MPI_COMM_NULL)
160  {
161  int i;
162  libmesh_call_mpi
163  (MPI_Comm_size(_communicator, &i));
164 
165  libmesh_assert_greater_equal (i, 0);
166  _size = cast_int<processor_id_type>(i);
167 
168  libmesh_call_mpi
169  (MPI_Comm_rank(_communicator, &i));
170 
171  libmesh_assert_greater_equal (i, 0);
172  _rank = cast_int<processor_id_type>(i);
173  }
174  else
175  {
176  _rank = 0;
177  _size = 1;
178  }
179 #endif
181 }
182 
183 
187 #ifdef LIBMESH_HAVE_MPI
189 {
190  if (this->size() > 1)
191  {
192  LOG_SCOPE("barrier()", "Parallel");
193  libmesh_call_mpi(MPI_Barrier (this->get()));
194  }
195 }
196 #else
197 void Communicator::barrier () const {}
198 #endif
199 
200 
202 {
203  if (used_tag_values.count(tagvalue))
204  {
205  // Get the largest value in the used values, and pick one
206  // larger
207  tagvalue = used_tag_values.rbegin()->first+1;
208  libmesh_assert(!used_tag_values.count(tagvalue));
209  }
210  used_tag_values[tagvalue] = 1;
211 
212  // #ifndef NDEBUG
213  // // Make sure everyone called get_unique_tag and make sure
214  // // everyone got the same value
215  // int maxval = tagvalue;
216  // this->max(maxval);
217  // libmesh_assert_equal_to (tagvalue, maxval);
218  // #endif
219 
220  return MessageTag(tagvalue, this);
221 }
222 
223 
224 } // namespace Parallel
225 
226 } // namespace libMesh
processor_id_type size() const
Definition: communicator.h:175
MPI_Comm communicator
Definition: communicator.h:57
MessageTag get_unique_tag(int tagvalue) const
Definition: communicator.C:201
void send_mode(const SendMode sm)
Definition: communicator.h:205
void assign(const communicator &comm)
Definition: communicator.C:155
void duplicate(const Communicator &comm)
Definition: communicator.C:108
Communicator & operator=(const Communicator &)=delete
void dereference_unique_tag(int tagvalue) const
Definition: communicator.C:44
std::map< int, unsigned int > used_tag_values
Definition: communicator.h:196
void reference_unique_tag(int tagvalue) const
Definition: communicator.C:35
void split(int color, int key, Communicator &target) const
Definition: communicator.C:89