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 /*
2 // @HEADER
3 // ***********************************************************************
4 //
5 // Teuchos: Common Tools Package
6 // Copyright (2004) Sandia Corporation
7 //
8 // Under terms of Contract DE-AC04-94AL85000, there is a non-exclusive
9 // license for use of this work by or on behalf of the U.S. Government.
10 //
11 // Redistribution and use in source and binary forms, with or without
12 // modification, are permitted provided that the following conditions are
13 // met:
14 //
15 // 1. Redistributions of source code must retain the above copyright
16 // notice, this list of conditions and the following disclaimer.
17 //
18 // 2. Redistributions in binary form must reproduce the above copyright
19 // notice, this list of conditions and the following disclaimer in the
20 // documentation and/or other materials provided with the distribution.
21 //
22 // 3. Neither the name of the Corporation nor the names of the
23 // contributors may be used to endorse or promote products derived from
24 // this software without specific prior written permission.
25 //
26 // THIS SOFTWARE IS PROVIDED BY SANDIA CORPORATION "AS IS" AND ANY
27 // EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE
28 // IMPLIED WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR
29 // PURPOSE ARE DISCLAIMED. IN NO EVENT SHALL SANDIA CORPORATION OR THE
30 // CONTRIBUTORS BE LIABLE FOR ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL,
31 // EXEMPLARY, OR CONSEQUENTIAL DAMAGES (INCLUDING, BUT NOT LIMITED TO,
32 // PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES; LOSS OF USE, DATA, OR
33 // PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON ANY THEORY OF
34 // LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT (INCLUDING
35 // NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
36 // SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
37 //
38 // Questions? Contact Michael A. Heroux (maherou@sandia.gov)
39 //
40 // ***********************************************************************
41 // @HEADER
42 */
43 
44 #include "Teuchos_CommHelpers.hpp"
45 #include "Teuchos_DefaultComm.hpp"
50 #include "Teuchos_Version.hpp"
51 #include "Teuchos_ScalarTraits.hpp"
53 #include "Teuchos_TimeMonitor.hpp"
54 #include "Teuchos_as.hpp"
55 
56 
57 //
58 // Unit test for Teuchos::Comm
59 //
60 
61 template<typename Ordinal>
63  const Teuchos::Comm<Ordinal> &comm,
65  const bool result
66  )
67 {
68  *out << "\nChecking that the above test passed in all processes ...";
69  int thisResult = ( result ? 1 : 0 );
70  int sumResult = -1;
71  reduceAll(comm,Teuchos::REDUCE_SUM,Ordinal(1),&thisResult,&sumResult);
72  const bool passed = sumResult==size(comm);
73  if(passed)
74  *out << " passed\n";
75  else
76  *out << " (sumResult="<<sumResult<<"!=numProcs) failed\n";
77  return passed;
78 }
79 
80 
81 template<typename Ordinal, typename Packet>
82 bool testComm(
83  const Teuchos::Comm<Ordinal> &comm,
85  )
86 {
87  using Teuchos::RCP;
88  using Teuchos::rcp;
91  using Teuchos::OSTab;
92  using Teuchos::dyn_cast;
93  using Teuchos::as;
94 
97 
98  OSTab tab(out);
99 
100  bool success = true, result;
101 
102  *out
103  << "\n***"
104  << "\n*** testComm<"<<OT::name()<<","<<ST::name()<<">(...)"
105  << "\n***\n";
106 
107  *out << "\nTesting Comm = " << comm.description() << "\n";
108 
109  const int procRank = rank(comm);
110  const int numProcs = size(comm);
111 
112  *out
113  << "\nnumProcs = size(comm) = " << numProcs << "\n"
114  << "\nprocRank = rank(comm) = " << procRank << "\n";
115 
116  const Ordinal count = numProcs*2;
117 
118  Teuchos::Array<Packet> sendBuff(count), recvBuff(count), recvBuff2(count);
119  for( int i = 0; i < count; ++i )
120  sendBuff[i] = Packet(procRank+1)*Packet(i);
121 
122  //
123  // send/receive
124  //
125 
126  if(numProcs > 1) {
127 
128 #ifdef TEUCHOS_MPI_COMM_DUMP
129  Teuchos::MpiComm<Ordinal>::show_dump = true;
130 #endif
131 
132  if(procRank==numProcs-1) {
133  *out << "\nSending data from p="<<procRank<<" to the root process (see p=0 output!) ...\n";
134  send(comm,count,&sendBuff[0],0);
135  }
136 
137  if(procRank==0) {
138  *out << "\nReceiving data specifically from p="<<numProcs-1<<" ...\n";
139  std::fill_n(&recvBuff[0],count,Packet(0));
140  const int sourceRank = receive(comm,numProcs-1,count,&recvBuff[0]);
141  result = sourceRank ==numProcs-1;
142  *out
143  << "\nChecking that sourceRank="<<sourceRank<<" == numProcs-1="<<(numProcs-1)
144  << " : " << (result ? "passed" : "falied" ) << "\n";
145  *out << "\nChecking that recvBuffer[] == numProcs * sendBuffer[] ...";
146  result = true;
147  for( int i = 0; i < count; ++i ) {
148  const Packet expected = Packet(numProcs)*sendBuff[i];
149  if( recvBuff[i] != expected ) {
150  result = false;
151  *out
152  << "\n recvBuffer["<<i<<"]="<<recvBuff[i]
153  << " == numProcs*sendBuffer["<<i<<"]="<<expected<<" : failed";
154  }
155  }
156  if(result) {
157  *out << " passed\n";
158  }
159  else {
160  *out << "\n";
161  success = false;
162  }
163  }
164 
165 #ifdef TEUCHOS_MPI_COMM_DUMP
166  Teuchos::MpiComm<Ordinal>::show_dump = false;
167 #endif
168 
169  }
170 
171 
172  //
173  // broadcast/reduceAll(sum)
174  //
175 
176  if(procRank==0) {
177  std::copy(&sendBuff[0],&sendBuff[0]+count,&recvBuff[0]);
178  *out << "\nSending broadcast of data from sendBuff[] in root process to recvBuff[] in each process ...\n";
179  }
180  else {
181  std::fill_n(&recvBuff[0],count,Packet(0));
182  *out << "\nReceiving broadcast of data from root process into recvBuff[] ...\n";
183  }
184 
185  broadcast(comm,0,count,&recvBuff[0]);
186 
187  *out << "\nSumming broadcasted data recvBuff[] over all processes into recvBuff2[] ...\n";
188 
189  reduceAll(comm,Teuchos::REDUCE_SUM,count,&recvBuff[0],&recvBuff2[0]);
190 
191  *out << "\nChecking that recvBuff2[i] == numProcs * i ...";
192  result = true;
193  for( int i = 0; i < count; ++i ) {
194  const Packet expected = Packet(numProcs)*Packet(i);
195  //*out << "\nexpected["<<i<<"]=numProcs*i="<<Packet(numProcs)<<"*"<<Packet(i)<<"="<<expected<<"\n";
196  if( recvBuff2[i] != expected ) {
197  result = false;
198  *out
199  << "\n recvBuffer2["<<i<<"]="<<recvBuff2[i]
200  << " == numProcs*"<<i<<"="<<expected<<" : failed";
201  }
202  }
203  if(result) {
204  *out << " passed\n";
205  }
206  else {
207  *out << "\n";
208  success = false;
209  }
210 
211  result = checkSumResult(comm,out,result);
212  if(!result) success = false;
213 
214  //
215  // reduceAll(min)
216  //
217 
218  if( ST::isComparable ) {
219 
220  *out << "\nTaking min of sendBuff[] and putting it in recvBuff[] ...\n";
221 
222  reduceAll(comm,Teuchos::REDUCE_MIN,count,&sendBuff[0],&recvBuff[0]);
223 
224  *out << "\nChecking that recvBuff[i] == i ...";
225  result = true;
226  for( int i = 0; i < count; ++i ) {
227  const Packet expected = Packet(i);
228  //*out << "\nexpected["<<i<<"]=numProcs*i="<<Packet(numProcs)<<"*"<<Packet(i)<<"="<<expected<<"\n";
229  if( recvBuff[i] != expected ) {
230  result = false;
231  *out
232  << "\n recvBuffer["<<i<<"]="<<recvBuff[i]
233  << " == "<<i<<"="<<expected<<" : failed";
234  }
235  }
236  if(result) {
237  *out << " passed\n";
238  }
239  else {
240  *out << "\n";
241  success = false;
242  }
243 
244  result = checkSumResult(comm,out,result);
245  if(!result) success = false;
246 
247  }
248 
249  //
250  // reduceAll(max)
251  //
252 
253  if( ST::isComparable ) {
254 
255  *out << "\nTaking max of sendBuff[] and putting it in recvBuff[] ...\n";
256 
257  reduceAll(comm,Teuchos::REDUCE_MAX,count,&sendBuff[0],&recvBuff[0]);
258 
259  *out << "\nChecking that recvBuff[i] == numProcs*i ...";
260  result = true;
261  for( int i = 0; i < count; ++i ) {
262  const Packet expected = Packet(numProcs)*Packet(i);
263  //*out << "\nexpected["<<i<<"]=numProcs*i="<<Packet(numProcs)<<"*"<<Packet(i)<<"="<<expected<<"\n";
264  if( recvBuff[i] != expected ) {
265  result = false;
266  *out
267  << "\n recvBuffer["<<i<<"]="<<recvBuff[i]
268  << " == numProcs*"<<i<<"="<<expected<<" : failed";
269  }
270  }
271  if(result) {
272  *out << " passed\n";
273  }
274  else {
275  *out << "\n";
276  success = false;
277  }
278 
279  result = checkSumResult(comm,out,result);
280  if(!result) success = false;
281 
282  }
283 
284  //
285  // gatherAll
286  //
287 
288  *out << "\nGathering all data from sendBuff[] in each process to all processes to allRecvBuff ...\n";
289 
291  allRecvBuff(count*numProcs);
292 
293  gatherAll(comm,count,&sendBuff[0],Ordinal(allRecvBuff.size()),&allRecvBuff[0]);
294 
295  *out << "\nChecking that allRecvBuff[count*k+i] == (k+1) * i ...";
296  result = true;
297  for( int k = 0; k < numProcs; ++k ) {
298  for( int i = 0; i < count; ++i ) {
299  const Packet expected = Packet(k+1)*Packet(i);
300  if( allRecvBuff[count*k+i] != expected ) {
301  result = false;
302  *out
303  << "\n allRecvBuff["<<count<<"*"<<k<<"+"<<i<<"]="<<allRecvBuff[count*k+i]
304  << " == (k+1)*i="<<expected<<" : failed";
305  }
306  }
307  }
308  if(result) {
309  *out << " passed\n";
310  }
311  else {
312  *out << "\n";
313  success = false;
314  }
315 
316  result = checkSumResult(comm,out,result);
317  if(!result) success = false;
318 
319  //
320  // scan
321  //
322 
323  *out << "\nPerforming a scan sum of sendBuff[] into recvBuff[] ...\n";
324 
325  std::fill_n(&recvBuff[0],count,Packet(0));
326 
327  scan(comm,Teuchos::REDUCE_SUM,count,&sendBuff[0],&recvBuff[0]);
328 
329  *out << "\nChecking that recvBuff[i] == sum(k+1,k=0...procRank) * i ...";
330  result = true;
331  int sumProcRank = 0;
332  for( int k = 0; k <= procRank; ++k ) sumProcRank += (k+1);
333  for( int i = 0; i < count; ++i ) {
334  const Packet expected = Packet(sumProcRank)*Packet(i);
335  //*out << "\nexpected["<<i<<"]=sum(k+1,k=0...procRank)*i="<<Packet(sumProcRank)<<"*"<<Packet(i)<<"="<<expected<<"\n";
336  if( recvBuff[i] != expected ) {
337  result = false;
338  *out
339  << "\n recvBuffer["<<i<<"]="<<recvBuff[i]
340  << " == sum(k+1,k=0...procRank)*"<<i<<"="<<expected<<" : failed";
341  }
342  }
343  if(result) {
344  *out << " passed\n";
345  }
346  else {
347  *out << "\n";
348  success = false;
349  }
350 
351  result = checkSumResult(comm,out,result);
352  if(!result) success = false;
353 
354  //
355  // The End!
356  //
357 
358  if(success)
359  *out << "\nCongratulations, all tests for this Comm check out!\n";
360  else
361  *out << "\nOh no, at least one of the tests for this Comm failed!\n";
362 
363  return success;
364 
365 }
366 
367 template<typename Ordinal>
370  )
371 {
372 
373  bool success = true, result;
374 
375  using Teuchos::RCP;
376  using Teuchos::rcp;
377  using Teuchos::FancyOStream;
379  using Teuchos::OSTab;
380 
382 
383  OSTab tab(out);
384 
385  RCP<const Teuchos::Comm<Ordinal> >
387 
388 #ifdef HAVE_MPI
389 
390  // Test that the DefaultComm is really a DefaultMpiComm.
391  RCP<const Teuchos::MpiComm<Ordinal> >
392  mpiComm = Teuchos::rcp_dynamic_cast<const Teuchos::MpiComm<Ordinal> >( comm, false );
393 
394  if (mpiComm == Teuchos::null) {
395  success = false;
396  *out << "\n*** FAILED to cast the Teuchos::DefaultComm<"<< OT::name() << "> to a Teuchos::MpiComm<" << OT::name() << ">!\n";
397  }
398  else {
399  *out
400  << "\n***"
401  << "\n*** Successfully casted the Teuchos::DefaultComm<"<< OT::name() << "> to a Teuchos::MpiComm<" << OT::name() << ">!"
402  << "\n***\n";
403 
404  // Now get the raw pointer to the MPI_Comm object
405  RCP<const Teuchos::OpaqueWrapper<MPI_Comm> >
406  rawMpiComm = mpiComm->getRawMpiComm();
407 
408  if (static_cast<MPI_Comm>(*rawMpiComm) == 0) {
409  success = false;
410  *out << "\n*** FAILED to get the raw MPI_Comm pointer from the Teuchos::MpiComm<" << OT::name() << ">!\n";
411  }
412  else {
413  *out
414  << "\n***"
415  << "\n*** Successfully got the raw MPI_Comm pointer from the Teuchos::MpiComm<" << OT::name() << ">!"
416  << "\n***\n";
417  }
418  }
419 
420 #endif
421 
422  *out
423  << "\n***"
424  << "\n*** Created a Comm of type " << comm->description() << " for testing"
425  << "\n***\n";
426 
427  *out << "\nOrdinal type = "<<OT::name()<<" with an extent of "<<sizeof(Ordinal)<<" bytes\n";
428 
429  if( comm->getSize() <= 4 ) {
430  result = testComm<Ordinal,char>(*comm,out);
431  if(!result) success = false;
432  }
433 
434  result = testComm<Ordinal,int>(*comm,out);
435  if(!result) success = false;
436 
437  result = testComm<Ordinal,size_t>(*comm,out);
438  if(!result) success = false;
439 
440  result = testComm<Ordinal,float>(*comm,out);
441  if(!result) success = false;
442 
443  result = testComm<Ordinal,double>(*comm,out);
444  if(!result) success = false;
445 
446 #ifdef HAVE_TEUCHOS_COMPLEX
447 
448 # ifdef HAVE_TEUCHOS_FLOAT
449 
450  result = testComm<Ordinal,std::complex<float> >(*comm,out);
451  if(!result) success = false;
452 
453 # endif // HAVE_TEUCHOS_FLOAT
454 
455  result = testComm<Ordinal,std::complex<double> >(*comm,out);
456  if(!result) success = false;
457 
458 #endif // HAVE_TEUCHOS_COMPLEX
459 
460  return success;
461 
462 }
463 
464 //
465 // Main driver program
466 //
467 
468 int main(int argc, char* argv[])
469 {
470 
471  using Teuchos::RCP;
472  using Teuchos::rcp;
473  using Teuchos::FancyOStream;
475  using Teuchos::OSTab;
477 
478  bool success = true, result;
479 
480  Teuchos::GlobalMPISession mpiSession(&argc,&argv);
481 
482  try {
483 
484  CommandLineProcessor clp;
485  clp.throwExceptions(false);
486  clp.addOutputSetupOptions(true);
487 
488  bool showTimers = true;
489 
490  clp.setOption( "show-timers", "no-show-timers", &showTimers, "Determine if timers are shown or not" );
491 
492  CommandLineProcessor::EParseCommandLineReturn
493  parse_return = clp.parse(argc,argv);
494  if( parse_return != CommandLineProcessor::PARSE_SUCCESSFUL )
495  return parse_return;
496 
497  RCP<FancyOStream>
498  out = VerboseObjectBase::getDefaultOStream();
499 
500  *out << std::endl << Teuchos::Teuchos_Version() << std::endl << std::endl;
501 
502  result = masterTestComm<short int>(out);
503  if(!result) success = false;
504 
505  result = masterTestComm<int>(out);
506  if(!result) success = false;
507 
508  result = masterTestComm<long int>(out);
509  if(!result) success = false;
510 
511  if(showTimers) {
513  *out<<"\n"
514  ,out->getOutputToRootOnly() < 0 // Show local time or not
515  );
516  }
517 
518  if(success)
519  *out << "\nEnd Result: TEST PASSED\n";
520 
521  }
522  TEUCHOS_STANDARD_CATCH_STATEMENTS(true,std::cerr,success);
523 
524  return ( success ? 0 : 1 );
525 
526 }
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:82
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:62
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:368
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...