27 template<
typename Ordinal>
34 *out <<
"\nChecking that the above test passed in all processes ...";
35 int thisResult = ( result ? 1 : 0 );
38 const bool passed = sumResult==
size(comm);
42 *out <<
" (sumResult="<<sumResult<<
"!=numProcs) failed\n";
47 template<
typename Ordinal,
typename Packet>
66 bool success =
true, result;
70 <<
"\n*** testComm<"<<OT::name()<<
","<<ST::name()<<
">(...)"
73 *out <<
"\nTesting Comm = " << comm.
description() <<
"\n";
75 const int procRank =
rank(comm);
76 const int numProcs =
size(comm);
79 <<
"\nnumProcs = size(comm) = " << numProcs <<
"\n"
80 <<
"\nprocRank = rank(comm) = " << procRank <<
"\n";
82 const Ordinal count = numProcs*2;
85 for(
int i = 0; i < count; ++i )
86 sendBuff[i] = Packet(procRank+1)*Packet(i);
94 #ifdef TEUCHOS_MPI_COMM_DUMP
95 Teuchos::MpiComm<Ordinal>::show_dump =
true;
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);
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;
109 <<
"\nChecking that sourceRank="<<sourceRank<<
" == numProcs-1="<<(numProcs-1)
110 <<
" : " << (result ?
"passed" :
"falied" ) <<
"\n";
111 *out <<
"\nChecking that recvBuffer[] == numProcs * sendBuffer[] ...";
113 for(
int i = 0; i < count; ++i ) {
114 const Packet expected = Packet(numProcs)*sendBuff[i];
115 if( recvBuff[i] != expected ) {
118 <<
"\n recvBuffer["<<i<<
"]="<<recvBuff[i]
119 <<
" == numProcs*sendBuffer["<<i<<
"]="<<expected<<
" : failed";
131 #ifdef TEUCHOS_MPI_COMM_DUMP
132 Teuchos::MpiComm<Ordinal>::show_dump =
false;
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";
147 std::fill_n(&recvBuff[0],count,Packet(0));
148 *out <<
"\nReceiving broadcast of data from root process into recvBuff[] ...\n";
153 *out <<
"\nSumming broadcasted data recvBuff[] over all processes into recvBuff2[] ...\n";
157 *out <<
"\nChecking that recvBuff2[i] == numProcs * i ...";
159 for(
int i = 0; i < count; ++i ) {
160 const Packet expected = Packet(numProcs)*Packet(i);
162 if( recvBuff2[i] != expected ) {
165 <<
"\n recvBuffer2["<<i<<
"]="<<recvBuff2[i]
166 <<
" == numProcs*"<<i<<
"="<<expected<<
" : failed";
178 if(!result) success =
false;
184 if( ST::isComparable ) {
186 *out <<
"\nTaking min of sendBuff[] and putting it in recvBuff[] ...\n";
190 *out <<
"\nChecking that recvBuff[i] == i ...";
192 for(
int i = 0; i < count; ++i ) {
193 const Packet expected = Packet(i);
195 if( recvBuff[i] != expected ) {
198 <<
"\n recvBuffer["<<i<<
"]="<<recvBuff[i]
199 <<
" == "<<i<<
"="<<expected<<
" : failed";
211 if(!result) success =
false;
219 if( ST::isComparable ) {
221 *out <<
"\nTaking max of sendBuff[] and putting it in recvBuff[] ...\n";
225 *out <<
"\nChecking that recvBuff[i] == numProcs*i ...";
227 for(
int i = 0; i < count; ++i ) {
228 const Packet expected = Packet(numProcs)*Packet(i);
230 if( recvBuff[i] != expected ) {
233 <<
"\n recvBuffer["<<i<<
"]="<<recvBuff[i]
234 <<
" == numProcs*"<<i<<
"="<<expected<<
" : failed";
246 if(!result) success =
false;
254 *out <<
"\nGathering all data from sendBuff[] in each process to all processes to allRecvBuff ...\n";
257 allRecvBuff(count*numProcs);
259 gatherAll(comm,count,&sendBuff[0],Ordinal(allRecvBuff.
size()),&allRecvBuff[0]);
261 *out <<
"\nChecking that allRecvBuff[count*k+i] == (k+1) * i ...";
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 ) {
269 <<
"\n allRecvBuff["<<count<<
"*"<<k<<
"+"<<i<<
"]="<<allRecvBuff[count*k+i]
270 <<
" == (k+1)*i="<<expected<<
" : failed";
283 if(!result) success =
false;
289 *out <<
"\nPerforming a scan sum of sendBuff[] into recvBuff[] ...\n";
291 std::fill_n(&recvBuff[0],count,Packet(0));
295 *out <<
"\nChecking that recvBuff[i] == sum(k+1,k=0...procRank) * i ...";
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);
302 if( recvBuff[i] != expected ) {
305 <<
"\n recvBuffer["<<i<<
"]="<<recvBuff[i]
306 <<
" == sum(k+1,k=0...procRank)*"<<i<<
"="<<expected<<
" : failed";
318 if(!result) success =
false;
325 *out <<
"\nCongratulations, all tests for this Comm check out!\n";
327 *out <<
"\nOh no, at least one of the tests for this Comm failed!\n";
333 template<
typename Ordinal>
339 bool success =
true, result;
351 RCP<const Teuchos::Comm<Ordinal> >
357 RCP<const Teuchos::MpiComm<Ordinal> >
358 mpiComm = Teuchos::rcp_dynamic_cast<
const Teuchos::MpiComm<Ordinal> >( comm, false );
362 *out <<
"\n*** FAILED to cast the Teuchos::DefaultComm<"<< OT::name() <<
"> to a Teuchos::MpiComm<" << OT::name() <<
">!\n";
367 <<
"\n*** Successfully casted the Teuchos::DefaultComm<"<< OT::name() <<
"> to a Teuchos::MpiComm<" << OT::name() <<
">!"
371 RCP<const Teuchos::OpaqueWrapper<MPI_Comm> >
372 rawMpiComm = mpiComm->getRawMpiComm();
374 if (static_cast<MPI_Comm>(*rawMpiComm) == 0) {
376 *out <<
"\n*** FAILED to get the raw MPI_Comm pointer from the Teuchos::MpiComm<" << OT::name() <<
">!\n";
381 <<
"\n*** Successfully got the raw MPI_Comm pointer from the Teuchos::MpiComm<" << OT::name() <<
">!"
390 <<
"\n*** Created a Comm of type " << comm->description() <<
" for testing"
393 *out <<
"\nOrdinal type = "<<OT::name()<<
" with an extent of "<<
sizeof(Ordinal)<<
" bytes\n";
395 if( comm->getSize() <= 4 ) {
396 result = testComm<Ordinal,char>(*comm,out);
397 if(!result) success =
false;
400 result = testComm<Ordinal,int>(*comm,out);
401 if(!result) success =
false;
403 result = testComm<Ordinal,size_t>(*comm,out);
404 if(!result) success =
false;
406 result = testComm<Ordinal,float>(*comm,out);
407 if(!result) success =
false;
409 result = testComm<Ordinal,double>(*comm,out);
410 if(!result) success =
false;
412 #ifdef HAVE_TEUCHOS_COMPLEX
414 # ifdef HAVE_TEUCHOS_FLOAT
416 result = testComm<Ordinal,std::complex<float> >(*comm,out);
417 if(!result) success =
false;
419 # endif // HAVE_TEUCHOS_FLOAT
421 result = testComm<Ordinal,std::complex<double> >(*comm,out);
422 if(!result) success =
false;
424 #endif // HAVE_TEUCHOS_COMPLEX
434 int main(
int argc,
char* argv[])
444 bool success =
true, result;
450 CommandLineProcessor clp;
451 clp.throwExceptions(
false);
452 clp.addOutputSetupOptions(
true);
454 bool showTimers =
true;
456 clp.setOption(
"show-timers",
"no-show-timers", &showTimers,
"Determine if timers are shown or not" );
458 CommandLineProcessor::EParseCommandLineReturn
459 parse_return = clp.parse(argc,argv);
460 if( parse_return != CommandLineProcessor::PARSE_SUCCESSFUL )
464 out = VerboseObjectBase::getDefaultOStream();
468 result = masterTestComm<short int>(out);
469 if(!result) success =
false;
471 result = masterTestComm<int>(out);
472 if(!result) success =
false;
474 result = masterTestComm<long int>(out);
475 if(!result) success =
false;
480 ,out->getOutputToRootOnly() < 0
485 *out <<
"\nEnd Result: TEST PASSED\n";
490 return ( success ? 0 : 1 );
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)
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)
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<T&> 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.
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)
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...