Teuchos Package Browser (Single Doxygen Collection)  Version of the Day
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
Comm_test.cpp
Go to the documentation of this file.
1 // @HEADER
2 // *****************************************************************************
3 // Teuchos: Common Tools Package
4 //
5 // Copyright 2004 NTESS and the Teuchos contributors.
6 // SPDX-License-Identifier: BSD-3-Clause
7 // *****************************************************************************
8 // @HEADER
9 
10 #include "Teuchos_CommHelpers.hpp"
11 #include "Teuchos_DefaultComm.hpp"
16 #include "Teuchos_Version.hpp"
17 #include "Teuchos_ScalarTraits.hpp"
19 #include "Teuchos_TimeMonitor.hpp"
20 #include "Teuchos_as.hpp"
21 
22 
23 //
24 // Unit test for Teuchos::Comm
25 //
26 
27 template<typename Ordinal>
29  const Teuchos::Comm<Ordinal> &comm,
31  const bool result
32  )
33 {
34  *out << "\nChecking that the above test passed in all processes ...";
35  int thisResult = ( result ? 1 : 0 );
36  int sumResult = -1;
37  reduceAll(comm,Teuchos::REDUCE_SUM,Ordinal(1),&thisResult,&sumResult);
38  const bool passed = sumResult==size(comm);
39  if(passed)
40  *out << " passed\n";
41  else
42  *out << " (sumResult="<<sumResult<<"!=numProcs) failed\n";
43  return passed;
44 }
45 
46 
47 template<typename Ordinal, typename Packet>
48 bool testComm(
49  const Teuchos::Comm<Ordinal> &comm,
51  )
52 {
53  using Teuchos::RCP;
54  using Teuchos::rcp;
57  using Teuchos::OSTab;
58  using Teuchos::dyn_cast;
59  using Teuchos::as;
60 
63 
64  OSTab tab(out);
65 
66  bool success = true, result;
67 
68  *out
69  << "\n***"
70  << "\n*** testComm<"<<OT::name()<<","<<ST::name()<<">(...)"
71  << "\n***\n";
72 
73  *out << "\nTesting Comm = " << comm.description() << "\n";
74 
75  const int procRank = rank(comm);
76  const int numProcs = size(comm);
77 
78  *out
79  << "\nnumProcs = size(comm) = " << numProcs << "\n"
80  << "\nprocRank = rank(comm) = " << procRank << "\n";
81 
82  const Ordinal count = numProcs*2;
83 
84  Teuchos::Array<Packet> sendBuff(count), recvBuff(count), recvBuff2(count);
85  for( int i = 0; i < count; ++i )
86  sendBuff[i] = Packet(procRank+1)*Packet(i);
87 
88  //
89  // send/receive
90  //
91 
92  if(numProcs > 1) {
93 
94 #ifdef TEUCHOS_MPI_COMM_DUMP
95  Teuchos::MpiComm<Ordinal>::show_dump = true;
96 #endif
97 
98  if(procRank==numProcs-1) {
99  *out << "\nSending data from p="<<procRank<<" to the root process (see p=0 output!) ...\n";
100  send(comm,count,&sendBuff[0],0);
101  }
102 
103  if(procRank==0) {
104  *out << "\nReceiving data specifically from p="<<numProcs-1<<" ...\n";
105  std::fill_n(&recvBuff[0],count,Packet(0));
106  const int sourceRank = receive(comm,numProcs-1,count,&recvBuff[0]);
107  result = sourceRank ==numProcs-1;
108  *out
109  << "\nChecking that sourceRank="<<sourceRank<<" == numProcs-1="<<(numProcs-1)
110  << " : " << (result ? "passed" : "falied" ) << "\n";
111  *out << "\nChecking that recvBuffer[] == numProcs * sendBuffer[] ...";
112  result = true;
113  for( int i = 0; i < count; ++i ) {
114  const Packet expected = Packet(numProcs)*sendBuff[i];
115  if( recvBuff[i] != expected ) {
116  result = false;
117  *out
118  << "\n recvBuffer["<<i<<"]="<<recvBuff[i]
119  << " == numProcs*sendBuffer["<<i<<"]="<<expected<<" : failed";
120  }
121  }
122  if(result) {
123  *out << " passed\n";
124  }
125  else {
126  *out << "\n";
127  success = false;
128  }
129  }
130 
131 #ifdef TEUCHOS_MPI_COMM_DUMP
132  Teuchos::MpiComm<Ordinal>::show_dump = false;
133 #endif
134 
135  }
136 
137 
138  //
139  // broadcast/reduceAll(sum)
140  //
141 
142  if(procRank==0) {
143  std::copy(&sendBuff[0],&sendBuff[0]+count,&recvBuff[0]);
144  *out << "\nSending broadcast of data from sendBuff[] in root process to recvBuff[] in each process ...\n";
145  }
146  else {
147  std::fill_n(&recvBuff[0],count,Packet(0));
148  *out << "\nReceiving broadcast of data from root process into recvBuff[] ...\n";
149  }
150 
151  broadcast(comm,0,count,&recvBuff[0]);
152 
153  *out << "\nSumming broadcasted data recvBuff[] over all processes into recvBuff2[] ...\n";
154 
155  reduceAll(comm,Teuchos::REDUCE_SUM,count,&recvBuff[0],&recvBuff2[0]);
156 
157  *out << "\nChecking that recvBuff2[i] == numProcs * i ...";
158  result = true;
159  for( int i = 0; i < count; ++i ) {
160  const Packet expected = Packet(numProcs)*Packet(i);
161  //*out << "\nexpected["<<i<<"]=numProcs*i="<<Packet(numProcs)<<"*"<<Packet(i)<<"="<<expected<<"\n";
162  if( recvBuff2[i] != expected ) {
163  result = false;
164  *out
165  << "\n recvBuffer2["<<i<<"]="<<recvBuff2[i]
166  << " == numProcs*"<<i<<"="<<expected<<" : failed";
167  }
168  }
169  if(result) {
170  *out << " passed\n";
171  }
172  else {
173  *out << "\n";
174  success = false;
175  }
176 
177  result = checkSumResult(comm,out,result);
178  if(!result) success = false;
179 
180  //
181  // reduceAll(min)
182  //
183 
184  if( ST::isComparable ) {
185 
186  *out << "\nTaking min of sendBuff[] and putting it in recvBuff[] ...\n";
187 
188  reduceAll(comm,Teuchos::REDUCE_MIN,count,&sendBuff[0],&recvBuff[0]);
189 
190  *out << "\nChecking that recvBuff[i] == i ...";
191  result = true;
192  for( int i = 0; i < count; ++i ) {
193  const Packet expected = Packet(i);
194  //*out << "\nexpected["<<i<<"]=numProcs*i="<<Packet(numProcs)<<"*"<<Packet(i)<<"="<<expected<<"\n";
195  if( recvBuff[i] != expected ) {
196  result = false;
197  *out
198  << "\n recvBuffer["<<i<<"]="<<recvBuff[i]
199  << " == "<<i<<"="<<expected<<" : failed";
200  }
201  }
202  if(result) {
203  *out << " passed\n";
204  }
205  else {
206  *out << "\n";
207  success = false;
208  }
209 
210  result = checkSumResult(comm,out,result);
211  if(!result) success = false;
212 
213  }
214 
215  //
216  // reduceAll(max)
217  //
218 
219  if( ST::isComparable ) {
220 
221  *out << "\nTaking max of sendBuff[] and putting it in recvBuff[] ...\n";
222 
223  reduceAll(comm,Teuchos::REDUCE_MAX,count,&sendBuff[0],&recvBuff[0]);
224 
225  *out << "\nChecking that recvBuff[i] == numProcs*i ...";
226  result = true;
227  for( int i = 0; i < count; ++i ) {
228  const Packet expected = Packet(numProcs)*Packet(i);
229  //*out << "\nexpected["<<i<<"]=numProcs*i="<<Packet(numProcs)<<"*"<<Packet(i)<<"="<<expected<<"\n";
230  if( recvBuff[i] != expected ) {
231  result = false;
232  *out
233  << "\n recvBuffer["<<i<<"]="<<recvBuff[i]
234  << " == numProcs*"<<i<<"="<<expected<<" : failed";
235  }
236  }
237  if(result) {
238  *out << " passed\n";
239  }
240  else {
241  *out << "\n";
242  success = false;
243  }
244 
245  result = checkSumResult(comm,out,result);
246  if(!result) success = false;
247 
248  }
249 
250  //
251  // gatherAll
252  //
253 
254  *out << "\nGathering all data from sendBuff[] in each process to all processes to allRecvBuff ...\n";
255 
257  allRecvBuff(count*numProcs);
258 
259  gatherAll(comm,count,&sendBuff[0],Ordinal(allRecvBuff.size()),&allRecvBuff[0]);
260 
261  *out << "\nChecking that allRecvBuff[count*k+i] == (k+1) * i ...";
262  result = true;
263  for( int k = 0; k < numProcs; ++k ) {
264  for( int i = 0; i < count; ++i ) {
265  const Packet expected = Packet(k+1)*Packet(i);
266  if( allRecvBuff[count*k+i] != expected ) {
267  result = false;
268  *out
269  << "\n allRecvBuff["<<count<<"*"<<k<<"+"<<i<<"]="<<allRecvBuff[count*k+i]
270  << " == (k+1)*i="<<expected<<" : failed";
271  }
272  }
273  }
274  if(result) {
275  *out << " passed\n";
276  }
277  else {
278  *out << "\n";
279  success = false;
280  }
281 
282  result = checkSumResult(comm,out,result);
283  if(!result) success = false;
284 
285  //
286  // scan
287  //
288 
289  *out << "\nPerforming a scan sum of sendBuff[] into recvBuff[] ...\n";
290 
291  std::fill_n(&recvBuff[0],count,Packet(0));
292 
293  scan(comm,Teuchos::REDUCE_SUM,count,&sendBuff[0],&recvBuff[0]);
294 
295  *out << "\nChecking that recvBuff[i] == sum(k+1,k=0...procRank) * i ...";
296  result = true;
297  int sumProcRank = 0;
298  for( int k = 0; k <= procRank; ++k ) sumProcRank += (k+1);
299  for( int i = 0; i < count; ++i ) {
300  const Packet expected = Packet(sumProcRank)*Packet(i);
301  //*out << "\nexpected["<<i<<"]=sum(k+1,k=0...procRank)*i="<<Packet(sumProcRank)<<"*"<<Packet(i)<<"="<<expected<<"\n";
302  if( recvBuff[i] != expected ) {
303  result = false;
304  *out
305  << "\n recvBuffer["<<i<<"]="<<recvBuff[i]
306  << " == sum(k+1,k=0...procRank)*"<<i<<"="<<expected<<" : failed";
307  }
308  }
309  if(result) {
310  *out << " passed\n";
311  }
312  else {
313  *out << "\n";
314  success = false;
315  }
316 
317  result = checkSumResult(comm,out,result);
318  if(!result) success = false;
319 
320  //
321  // The End!
322  //
323 
324  if(success)
325  *out << "\nCongratulations, all tests for this Comm check out!\n";
326  else
327  *out << "\nOh no, at least one of the tests for this Comm failed!\n";
328 
329  return success;
330 
331 }
332 
333 template<typename Ordinal>
336  )
337 {
338 
339  bool success = true, result;
340 
341  using Teuchos::RCP;
342  using Teuchos::rcp;
343  using Teuchos::FancyOStream;
345  using Teuchos::OSTab;
346 
348 
349  OSTab tab(out);
350 
351  RCP<const Teuchos::Comm<Ordinal> >
353 
354 #ifdef HAVE_MPI
355 
356  // Test that the DefaultComm is really a DefaultMpiComm.
357  RCP<const Teuchos::MpiComm<Ordinal> >
358  mpiComm = Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<Ordinal> >( comm, false );
359 
360  if (mpiComm == Teuchos::null) {
361  success = false;
362  *out << "\n*** FAILED to cast the Teuchos::DefaultComm<"<< OT::name() << "> to a Teuchos::MpiComm<" << OT::name() << ">!\n";
363  }
364  else {
365  *out
366  << "\n***"
367  << "\n*** Successfully casted the Teuchos::DefaultComm<"<< OT::name() << "> to a Teuchos::MpiComm<" << OT::name() << ">!"
368  << "\n***\n";
369 
370  // Now get the raw pointer to the MPI_Comm object
371  RCP<const Teuchos::OpaqueWrapper<MPI_Comm> >
372  rawMpiComm = mpiComm->getRawMpiComm();
373 
374  if (static_cast<MPI_Comm>(*rawMpiComm) == 0) {
375  success = false;
376  *out << "\n*** FAILED to get the raw MPI_Comm pointer from the Teuchos::MpiComm<" << OT::name() << ">!\n";
377  }
378  else {
379  *out
380  << "\n***"
381  << "\n*** Successfully got the raw MPI_Comm pointer from the Teuchos::MpiComm<" << OT::name() << ">!"
382  << "\n***\n";
383  }
384  }
385 
386 #endif
387 
388  *out
389  << "\n***"
390  << "\n*** Created a Comm of type " << comm->description() << " for testing"
391  << "\n***\n";
392 
393  *out << "\nOrdinal type = "<<OT::name()<<" with an extent of "<<sizeof(Ordinal)<<" bytes\n";
394 
395  if( comm->getSize() <= 4 ) {
396  result = testComm<Ordinal,char>(*comm,out);
397  if(!result) success = false;
398  }
399 
400  result = testComm<Ordinal,int>(*comm,out);
401  if(!result) success = false;
402 
403  result = testComm<Ordinal,size_t>(*comm,out);
404  if(!result) success = false;
405 
406  result = testComm<Ordinal,float>(*comm,out);
407  if(!result) success = false;
408 
409  result = testComm<Ordinal,double>(*comm,out);
410  if(!result) success = false;
411 
412 #ifdef HAVE_TEUCHOS_COMPLEX
413 
414 # ifdef HAVE_TEUCHOS_FLOAT
415 
416  result = testComm<Ordinal,std::complex<float> >(*comm,out);
417  if(!result) success = false;
418 
419 # endif // HAVE_TEUCHOS_FLOAT
420 
421  result = testComm<Ordinal,std::complex<double> >(*comm,out);
422  if(!result) success = false;
423 
424 #endif // HAVE_TEUCHOS_COMPLEX
425 
426  return success;
427 
428 }
429 
430 //
431 // Main driver program
432 //
433 
434 int main(int argc, char* argv[])
435 {
436 
437  using Teuchos::RCP;
438  using Teuchos::rcp;
439  using Teuchos::FancyOStream;
441  using Teuchos::OSTab;
443 
444  bool success = true, result;
445 
446  Teuchos::GlobalMPISession mpiSession(&argc,&argv);
447 
448  try {
449 
450  CommandLineProcessor clp;
451  clp.throwExceptions(false);
452  clp.addOutputSetupOptions(true);
453 
454  bool showTimers = true;
455 
456  clp.setOption( "show-timers", "no-show-timers", &showTimers, "Determine if timers are shown or not" );
457 
458  CommandLineProcessor::EParseCommandLineReturn
459  parse_return = clp.parse(argc,argv);
460  if( parse_return != CommandLineProcessor::PARSE_SUCCESSFUL )
461  return parse_return;
462 
463  RCP<FancyOStream>
464  out = VerboseObjectBase::getDefaultOStream();
465 
466  *out << std::endl << Teuchos::Teuchos_Version() << std::endl << std::endl;
467 
468  result = masterTestComm<short int>(out);
469  if(!result) success = false;
470 
471  result = masterTestComm<int>(out);
472  if(!result) success = false;
473 
474  result = masterTestComm<long int>(out);
475  if(!result) success = false;
476 
477  if(showTimers) {
479  *out<<"\n"
480  ,out->getOutputToRootOnly() < 0 // Show local time or not
481  );
482  }
483 
484  if(success)
485  *out << "\nEnd Result: TEST PASSED\n";
486 
487  }
488  TEUCHOS_STANDARD_CATCH_STATEMENTS(true,std::cerr,success);
489 
490  return ( success ? 0 : 1 );
491 
492 }
void broadcast(const Comm< Ordinal > &comm, const int rootRank, const Ordinal count, Packet buffer[])
Broadcast array of objects that use value semantics.
bool testComm(const Teuchos::Comm< Ordinal > &comm, const Teuchos::RCP< Teuchos::FancyOStream > &out)
Definition: Comm_test.cpp:48
int rank(const Comm< Ordinal > &comm)
Get the process rank.
static Teuchos::RCP< const Comm< OrdinalType > > getComm()
Return the default global communicator.
Initialize, finalize, and query the global MPI session.
basic_OSTab< char > OSTab
void reduceAll(const Comm< Ordinal > &comm, const ValueTypeReductionOp< Ordinal, Packet > &reductOp, const Ordinal count, const Packet sendBuffer[], Packet globalReducts[])
Wrapper for MPI_Allreduce that takes a custom reduction operator.
This structure defines some basic traits for a scalar field type.
bool checkSumResult(const Teuchos::Comm< Ordinal > &comm, const Teuchos::RCP< Teuchos::FancyOStream > &out, const bool result)
Definition: Comm_test.cpp:28
TEUCHOS_DEPRECATED RCP< T > rcp(T *p, Dealloc_T dealloc, bool owns_mem)
Deprecated.
static void summarize(Ptr< const Comm< int > > comm, std::ostream &out=std::cout, const bool alwaysWriteLocal=false, const bool writeGlobalStats=true, const bool writeZeroTimers=true, const ECounterSetOp setOp=Intersection, const std::string &filter="", const bool ignoreZeroTimers=false)
Print summary statistics for all timers on the given communicator.
void send(const Comm< Ordinal > &comm, const Ordinal count, const Packet sendBuffer[], const int destRank)
Send objects that use values semantics to another process.
T_To & dyn_cast(T_From &from)
Dynamic casting utility function meant to replace dynamic_cast&lt;T&amp;&gt; by throwing a better documented er...
#define TEUCHOS_STANDARD_CATCH_STATEMENTS(VERBOSE, ERR_STREAM, SUCCESS_FLAG)
Simple macro that catches and reports standard exceptions and other exceptions.
virtual std::string description() const
Return a simple one-line description of this object.
Non-templated base class for objects that can print their activities to a stream. ...
This structure defines some basic traits for the ordinal field type.
std::string Teuchos_Version()
int size(const Comm< Ordinal > &comm)
Get the number of processes in the communicator.
TypeTo as(const TypeFrom &t)
Convert from one value type to another.
int main(int argc, char *argv[])
basic_FancyOStream< char > FancyOStream
Defines basic traits for the ordinal field type.
A MPI utilities class, providing methods for initializing, finalizing, and querying the global MPI se...
Basic command line parser for input from (argc,argv[])
int receive(const Comm< Ordinal > &comm, const int sourceRank, const Ordinal count, Packet recvBuffer[])
Receive objects that use values semantics from another process.
size_type size() const
Defines basic traits for the scalar field type.
Scope guard for Teuchos::Time, with MPI collective timer reporting.
Smart reference counting pointer class for automatic garbage collection.
bool masterTestComm(const Teuchos::RCP< Teuchos::FancyOStream > &out)
Definition: Comm_test.cpp:334
void gatherAll(const Comm< Ordinal > &comm, const Ordinal sendCount, const Packet sendBuffer[], const Ordinal recvCount, Packet recvBuffer[])
Gather array of objects that use value semantics from every process to every process.
RCP< basic_FancyOStream< CharT, Traits > > tab(const RCP< basic_FancyOStream< CharT, Traits > > &out, const int tabs=1, const std::basic_string< CharT, Traits > linePrefix="")
Create a tab for an RCP-wrapped basic_FancyOStream object to cause the indentation of all output auto...
Definition of Teuchos::as, for conversions between types.
Class that helps parse command line input arguments from (argc,argv[]) and set options.
void scan(const Comm< Ordinal > &comm, const ValueTypeReductionOp< Ordinal, Packet > &reductOp, const Ordinal count, const Packet sendBuffer[], Packet scanReducts[])
Scan/Reduce array of objects that use value semantics using a user-defined reduction operator...
Replacement for std::vector that is compatible with the Teuchos Memory Management classes...