queso-0.51.1
MpiComm.C
Go to the documentation of this file.
1 //-----------------------------------------------------------------------bl-
2 //--------------------------------------------------------------------------
3 //
4 // QUESO - a library to support the Quantification of Uncertainty
5 // for Estimation, Simulation and Optimization
6 //
7 // Copyright (C) 2008,2009,2010,2011,2012,2013 The PECOS Development Team
8 //
9 // This library is free software; you can redistribute it and/or
10 // modify it under the terms of the Version 2.1 GNU Lesser General
11 // Public License as published by the Free Software Foundation.
12 //
13 // This library is distributed in the hope that it will be useful,
14 // but WITHOUT ANY WARRANTY; without even the implied warranty of
15 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
16 // Lesser General Public License for more details.
17 //
18 // You should have received a copy of the GNU Lesser General Public
19 // License along with this library; if not, write to the Free Software
20 // Foundation, Inc. 51 Franklin Street, Fifth Floor,
21 // Boston, MA 02110-1301 USA
22 //
23 //-----------------------------------------------------------------------el-
24 
25 #include <queso/MpiComm.h>
26 #include <queso/Environment.h>
27 
28 namespace QUESO {
29 
30 // Default constructor ------------------------------
32  :
33  m_env( *(new EmptyEnvironment()) )
34 {
37  "MpiComm::constructor()",
38  "should not be called");
39 }
40 // QUESO MpiComm MPI Constructor ------------------
42  :
43  m_env (env),
44 #ifdef QUESO_HAS_TRILINOS
45  m_epetraMpiComm( new Epetra_MpiComm(inputRawComm) ),
46 #endif
47  m_rawComm (inputRawComm),
48  m_worldRank (-1),
49  m_myPid (-1),
50  m_numProc (-1)
51 {
52  int mpiRC = MPI_Comm_rank(inputRawComm,&m_worldRank);
53  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
55  "MpiComm::constructor()",
56  "failed MPI_Comm_rank() on full rank");
57 
58  mpiRC = MPI_Comm_rank(inputRawComm,&m_myPid);
59  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
61  "MpiComm::constructor()",
62  "failed MPI_Comm_rank() on inputRawComm");
63 
64  mpiRC = MPI_Comm_size(inputRawComm,&m_numProc);
65  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
67  "MpiComm::constructor()",
68  "failed MPI_Comm_size() on inputRawComm");
69 }
70 
71 // Copy constructor ---------------------------------
73  :
74  m_env (src.m_env)
75 #ifdef QUESO_HAS_TRILINOS
76  ,
77  m_epetraMpiComm(NULL)
78 #endif
79 {
80  this->copy(src);
81 }
82 
83 // Destructor ---------------------------------------
85 {
86 #ifdef QUESO_HAS_TRILINOS
87  delete m_epetraMpiComm;
88  m_epetraMpiComm = NULL;
89 #endif
90 }
91 
92 // --------------------------------------------------
93 // Set methodos -------------------------------------
94 MpiComm&
96 {
97  this->copy(rhs);
98  return *this;
99 }
100 
101 // Attribute access methods -------------------------
104 {
105 #ifdef QUESO_HAS_TRILINOS
106  return m_epetraMpiComm->Comm();
107 #endif
108  return m_rawComm;
109 }
110 // --------------------------------------------------
111 int
113 {
114 #ifdef QUESO_HAS_TRILINOS
115  return m_epetraMpiComm->MyPID();
116 #endif
117  return m_myPid;
118 }
119 // --------------------------------------------------
120 int
122 {
123 #ifdef QUESO_HAS_TRILINOS
124  return m_epetraMpiComm->NumProc();
125 #endif
126  return m_numProc;
127 }
128 // Methods overridden from Comm ---------------------
129 
130 void
131 MpiComm::Allreduce(void* sendbuf, void* recvbuf, int count, RawType_MPI_Datatype datatype, RawType_MPI_Op op, const char* whereMsg, const char* whatMsg) const
132 {
133  int mpiRC = MPI_Allreduce(sendbuf, recvbuf, count, datatype, op, m_rawComm);
134  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
135  m_worldRank,
136  whereMsg,
137  whatMsg);
138 
139  return;
140 }
141 //--------------------------------------------------
142 void
143 MpiComm::Barrier() const // const char* whereMsg, const char* whatMsg) const
144 {
145 #ifdef QUESO_HAS_TRILINOS
146  return m_epetraMpiComm->Barrier();
147 #endif
148  int mpiRC = MPI_Barrier(m_rawComm);
149  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
150  m_worldRank,
151  "MPIComm::Barrier()", // whereMsg,
152  "mpiRC indicates failure"); // whatMsg);
153  return;
154 }
155 //--------------------------------------------------
156 void
157 MpiComm::Bcast(void* buffer, int count, RawType_MPI_Datatype datatype, int root, const char* whereMsg, const char* whatMsg) const
158 {
159  int mpiRC = MPI_Bcast(buffer, count, datatype, root, m_rawComm);
160  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
161  m_worldRank,
162  whereMsg,
163  whatMsg);
164  return;
165 }
166 //--------------------------------------------------
167 void
169  void* sendbuf, int sendcnt, RawType_MPI_Datatype sendtype,
170  void* recvbuf, int recvcount, RawType_MPI_Datatype recvtype,
171  int root,
172  const char* whereMsg, const char* whatMsg) const
173 {
174  //int MPI_Gather (void *sendbuf, int sendcnt, MPI_Datatype sendtype,
175  // void *recvbuf, int recvcount, MPI_Datatype recvtype,
176  // int root, MPI_Comm comm )
177  int mpiRC = MPI_Gather(sendbuf, sendcnt, sendtype,
178  recvbuf, recvcount, recvtype,
179  root, m_rawComm);
180  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
181  m_worldRank,
182  whereMsg,
183  whatMsg);
184  return;
185 }
186 //--------------------------------------------------
187 void
189  void* sendbuf, int sendcnt, RawType_MPI_Datatype sendtype,
190  void* recvbuf, int* recvcnts, int* displs, RawType_MPI_Datatype recvtype,
191  int root,
192  const char* whereMsg, const char* whatMsg) const
193 {
194  //int MPI_Gatherv(void *sendbuf, int sendcnt, MPI_Datatype sendtype,
195  // void *recvbuf, int *recvcnts, int *displs, MPI_Datatype recvtype,
196  // int root, MPI_Comm comm )
197  int mpiRC = MPI_Gatherv(sendbuf, sendcnt, sendtype,
198  recvbuf, recvcnts, displs, recvtype,
199  root, m_rawComm);
200  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
201  m_worldRank,
202  whereMsg,
203  whatMsg);
204  return;
205 }
206 //--------------------------------------------------
207 void
209  void* buf, int count, RawType_MPI_Datatype datatype, int source, int tag, RawType_MPI_Status* status,
210  const char* whereMsg, const char* whatMsg) const
211 {
212  int mpiRC = MPI_Recv(buf, count, datatype, source, tag, m_rawComm, status);
213  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
214  m_worldRank,
215  whereMsg,
216  whatMsg);
217  return;
218 }
219 //--------------------------------------------------
220 void
222  void* buf, int count, RawType_MPI_Datatype datatype, int dest, int tag,
223  const char* whereMsg, const char* whatMsg) const
224 {
225  int mpiRC = MPI_Send(buf, count, datatype, dest, tag, m_rawComm);
226  UQ_FATAL_TEST_MACRO(mpiRC != MPI_SUCCESS,
227  m_worldRank,
228  whereMsg,
229  whatMsg);
230  return;
231 }
232 // Misc methods ------------------------------------
233 void
234 MpiComm::syncPrintDebugMsg(const char* msg, unsigned int msgVerbosity, unsigned int numUSecs) const
235 {
236  if (m_env.syncVerbosity() >= msgVerbosity) {
237  this->Barrier();
238  for (int i = 0; i < this->NumProc(); ++i) {
239  if (i == this->MyPID()) {
240  std::cout << msg
241  << ": fullRank " << m_env.fullRank()
242  << ", subEnvironment " << m_env.subId()
243  << ", subRank " << m_env.subRank()
244  << ", inter0Rank " << m_env.inter0Rank()
245  << std::endl;
246  }
247  usleep(numUSecs);
248  this->Barrier();
249  }
250  //if (this->fullRank() == 0) std::cout << "Sleeping " << numUSecs << " microseconds..."
251  // << std::endl;
252  //usleep(numUSecs);
253  this->Barrier();
254  }
255 
256  return;
257 }
258 // -------------------------------------------------
259 #ifdef QUESO_HAS_TRILINOS
260 const Epetra_MpiComm&
261 MpiComm::epetraMpiComm() const
262 {
263  return *m_epetraMpiComm;
264 }
265 #endif
266 
267 // Private methods----------------------------------
268 void
270 {
271 #ifdef QUESO_HAS_TRILINOS
272  delete m_epetraMpiComm;
273  m_epetraMpiComm = new Epetra_MpiComm(*src.m_epetraMpiComm);
274 #endif
275  m_rawComm = src.m_rawComm;
276  m_worldRank = src.m_worldRank;
277  m_myPid = src.m_myPid;
278  m_numProc = src.m_numProc;
279 
280  return;
281 }
282 // -------------------------------------------------
283 
284 } // End namespace QUESO
void Bcast(void *buffer, int count, RawType_MPI_Datatype datatype, int root, const char *whereMsg, const char *whatMsg) const
Broadcast values from the root process to the slave processes.
Definition: MpiComm.C:157
int m_worldRank
World rank.
Definition: MpiComm.h:216
The QUESO MPI Communicator Class.
Definition: MpiComm.h:75
int subRank() const
Access function for sub-rank.
Definition: Environment.C:263
RawType_MPI_Comm m_rawComm
Embedded wrapped opaque MPI_Comm object.
Definition: MpiComm.h:213
MpiComm()
Default Constructor.
Definition: MpiComm.C:31
unsigned int syncVerbosity() const
Access function to private attribute m_syncVerbosity.
Definition: Environment.C:446
int inter0Rank() const
Returns the process inter0 rank.
Definition: Environment.C:289
const BaseEnvironment & m_env
Definition: MpiComm.h:206
int NumProc() const
Returns total number of processes.
Definition: MpiComm.C:121
int m_myPid
Process ID of this process.
Definition: MpiComm.h:219
MPI_Datatype RawType_MPI_Datatype
Definition: MpiComm.h:40
This class sets up the environment underlying the use of the QUESO library by an executable.
Definition: Environment.h:396
MPI_Status RawType_MPI_Status
Definition: MpiComm.h:42
void Allreduce(void *sendbuf, void *recvbuf, int count, RawType_MPI_Datatype datatype, RawType_MPI_Op op, const char *whereMsg, const char *whatMsg) const
Combines values from all processes and distributes the result back to all processes.
Definition: MpiComm.C:131
const int UQ_UNAVAILABLE_RANK
Definition: Defines.h:74
This (virtual) class sets up the environment underlying the use of the QUESO library by an executable...
Definition: Environment.h:187
void Send(void *buf, int count, RawType_MPI_Datatype datatype, int dest, int tag, const char *whereMsg, const char *whatMsg) const
Possibly blocking send of data from this process to another process.
Definition: MpiComm.C:221
void Gather(void *sendbuf, int sendcnt, RawType_MPI_Datatype sendtype, void *recvbuf, int recvcount, RawType_MPI_Datatype recvtype, int root, const char *whereMsg, const char *whatMsg) const
Gather values from each process to collect on all processes.
Definition: MpiComm.C:168
void Barrier() const
Pause every process in *this communicator until all the processes reach this point.
Definition: MpiComm.C:143
MpiComm & operator=(const MpiComm &rhs)
Assignment operator.
Definition: MpiComm.C:95
MPI_Comm RawType_MPI_Comm
Definition: MpiComm.h:38
void Gatherv(void *sendbuf, int sendcnt, RawType_MPI_Datatype sendtype, void *recvbuf, int *recvcnts, int *displs, RawType_MPI_Datatype recvtype, int root, const char *whereMsg, const char *whatMsg) const
Gathers into specified locations from all processes in a group.
Definition: MpiComm.C:188
int fullRank() const
Returns the process full rank.
Definition: Environment.C:241
void syncPrintDebugMsg(const char *msg, unsigned int msgVerbosity, unsigned int numUSecs) const
Synchronizes all the processes and print debug message.
Definition: MpiComm.C:234
int MyPID() const
Return my process ID.
Definition: MpiComm.C:112
~MpiComm()
Destructor.
Definition: MpiComm.C:84
void copy(const MpiComm &src)
Copies from an existing MpiComm instance.
Definition: MpiComm.C:269
void Recv(void *buf, int count, RawType_MPI_Datatype datatype, int source, int tag, RawType_MPI_Status *status, const char *whereMsg, const char *whatMsg) const
Blocking receive of data from this process to another process.
Definition: MpiComm.C:208
#define UQ_FATAL_TEST_MACRO(test, givenRank, where, what)
Definition: Defines.h:223
unsigned int subId() const
Access function to the number of each sub-environment Id: m_subId.
Definition: Environment.C:329
MPI_Op RawType_MPI_Op
Definition: MpiComm.h:41
RawType_MPI_Comm Comm() const
Extract MPI Communicator from a MpiComm object.
Definition: MpiComm.C:103

Generated on Thu Apr 23 2015 19:26:15 for queso-0.51.1 by  doxygen 1.8.5