61 template<
typename Ordinal>
68 *out <<
"\nChecking that the above test passed in all processes ...";
69 int thisResult = ( result ? 1 : 0 );
72 const bool passed = sumResult==
size(comm);
76 *out <<
" (sumResult="<<sumResult<<
"!=numProcs) failed\n";
81 template<
typename Ordinal,
typename Packet>
100 bool success =
true, result;
104 <<
"\n*** testComm<"<<OT::name()<<
","<<ST::name()<<
">(...)"
107 *out <<
"\nTesting Comm = " << comm.
description() <<
"\n";
109 const int procRank =
rank(comm);
110 const int numProcs =
size(comm);
113 <<
"\nnumProcs = size(comm) = " << numProcs <<
"\n"
114 <<
"\nprocRank = rank(comm) = " << procRank <<
"\n";
116 const Ordinal count = numProcs*2;
119 for(
int i = 0; i < count; ++i )
120 sendBuff[i] = Packet(procRank+1)*Packet(i);
128 #ifdef TEUCHOS_MPI_COMM_DUMP
129 Teuchos::MpiComm<Ordinal>::show_dump =
true;
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);
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;
143 <<
"\nChecking that sourceRank="<<sourceRank<<
" == numProcs-1="<<(numProcs-1)
144 <<
" : " << (result ?
"passed" :
"falied" ) <<
"\n";
145 *out <<
"\nChecking that recvBuffer[] == numProcs * sendBuffer[] ...";
147 for(
int i = 0; i < count; ++i ) {
148 const Packet expected = Packet(numProcs)*sendBuff[i];
149 if( recvBuff[i] != expected ) {
152 <<
"\n recvBuffer["<<i<<
"]="<<recvBuff[i]
153 <<
" == numProcs*sendBuffer["<<i<<
"]="<<expected<<
" : failed";
165 #ifdef TEUCHOS_MPI_COMM_DUMP
166 Teuchos::MpiComm<Ordinal>::show_dump =
false;
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";
181 std::fill_n(&recvBuff[0],count,Packet(0));
182 *out <<
"\nReceiving broadcast of data from root process into recvBuff[] ...\n";
187 *out <<
"\nSumming broadcasted data recvBuff[] over all processes into recvBuff2[] ...\n";
191 *out <<
"\nChecking that recvBuff2[i] == numProcs * i ...";
193 for(
int i = 0; i < count; ++i ) {
194 const Packet expected = Packet(numProcs)*Packet(i);
196 if( recvBuff2[i] != expected ) {
199 <<
"\n recvBuffer2["<<i<<
"]="<<recvBuff2[i]
200 <<
" == numProcs*"<<i<<
"="<<expected<<
" : failed";
212 if(!result) success =
false;
218 if( ST::isComparable ) {
220 *out <<
"\nTaking min of sendBuff[] and putting it in recvBuff[] ...\n";
224 *out <<
"\nChecking that recvBuff[i] == i ...";
226 for(
int i = 0; i < count; ++i ) {
227 const Packet expected = Packet(i);
229 if( recvBuff[i] != expected ) {
232 <<
"\n recvBuffer["<<i<<
"]="<<recvBuff[i]
233 <<
" == "<<i<<
"="<<expected<<
" : failed";
245 if(!result) success =
false;
253 if( ST::isComparable ) {
255 *out <<
"\nTaking max of sendBuff[] and putting it in recvBuff[] ...\n";
259 *out <<
"\nChecking that recvBuff[i] == numProcs*i ...";
261 for(
int i = 0; i < count; ++i ) {
262 const Packet expected = Packet(numProcs)*Packet(i);
264 if( recvBuff[i] != expected ) {
267 <<
"\n recvBuffer["<<i<<
"]="<<recvBuff[i]
268 <<
" == numProcs*"<<i<<
"="<<expected<<
" : failed";
280 if(!result) success =
false;
288 *out <<
"\nGathering all data from sendBuff[] in each process to all processes to allRecvBuff ...\n";
291 allRecvBuff(count*numProcs);
293 gatherAll(comm,count,&sendBuff[0],Ordinal(allRecvBuff.
size()),&allRecvBuff[0]);
295 *out <<
"\nChecking that allRecvBuff[count*k+i] == (k+1) * i ...";
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 ) {
303 <<
"\n allRecvBuff["<<count<<
"*"<<k<<
"+"<<i<<
"]="<<allRecvBuff[count*k+i]
304 <<
" == (k+1)*i="<<expected<<
" : failed";
317 if(!result) success =
false;
323 *out <<
"\nPerforming a scan sum of sendBuff[] into recvBuff[] ...\n";
325 std::fill_n(&recvBuff[0],count,Packet(0));
329 *out <<
"\nChecking that recvBuff[i] == sum(k+1,k=0...procRank) * i ...";
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);
336 if( recvBuff[i] != expected ) {
339 <<
"\n recvBuffer["<<i<<
"]="<<recvBuff[i]
340 <<
" == sum(k+1,k=0...procRank)*"<<i<<
"="<<expected<<
" : failed";
352 if(!result) success =
false;
359 *out <<
"\nCongratulations, all tests for this Comm check out!\n";
361 *out <<
"\nOh no, at least one of the tests for this Comm failed!\n";
367 template<
typename Ordinal>
373 bool success =
true, result;
385 RCP<const Teuchos::Comm<Ordinal> >
391 RCP<const Teuchos::MpiComm<Ordinal> >
392 mpiComm = Teuchos::rcp_dynamic_cast<
const Teuchos::MpiComm<Ordinal> >( comm, false );
396 *out <<
"\n*** FAILED to cast the Teuchos::DefaultComm<"<< OT::name() <<
"> to a Teuchos::MpiComm<" << OT::name() <<
">!\n";
401 <<
"\n*** Successfully casted the Teuchos::DefaultComm<"<< OT::name() <<
"> to a Teuchos::MpiComm<" << OT::name() <<
">!"
405 RCP<const Teuchos::OpaqueWrapper<MPI_Comm> >
406 rawMpiComm = mpiComm->getRawMpiComm();
408 if (static_cast<MPI_Comm>(*rawMpiComm) == 0) {
410 *out <<
"\n*** FAILED to get the raw MPI_Comm pointer from the Teuchos::MpiComm<" << OT::name() <<
">!\n";
415 <<
"\n*** Successfully got the raw MPI_Comm pointer from the Teuchos::MpiComm<" << OT::name() <<
">!"
424 <<
"\n*** Created a Comm of type " << comm->description() <<
" for testing"
427 *out <<
"\nOrdinal type = "<<OT::name()<<
" with an extent of "<<
sizeof(Ordinal)<<
" bytes\n";
429 if( comm->getSize() <= 4 ) {
430 result = testComm<Ordinal,char>(*comm,out);
431 if(!result) success =
false;
434 result = testComm<Ordinal,int>(*comm,out);
435 if(!result) success =
false;
437 result = testComm<Ordinal,size_t>(*comm,out);
438 if(!result) success =
false;
440 result = testComm<Ordinal,float>(*comm,out);
441 if(!result) success =
false;
443 result = testComm<Ordinal,double>(*comm,out);
444 if(!result) success =
false;
446 #ifdef HAVE_TEUCHOS_COMPLEX
448 # ifdef HAVE_TEUCHOS_FLOAT
450 result = testComm<Ordinal,std::complex<float> >(*comm,out);
451 if(!result) success =
false;
453 # endif // HAVE_TEUCHOS_FLOAT
455 result = testComm<Ordinal,std::complex<double> >(*comm,out);
456 if(!result) success =
false;
458 #endif // HAVE_TEUCHOS_COMPLEX
468 int main(
int argc,
char* argv[])
478 bool success =
true, result;
484 CommandLineProcessor clp;
485 clp.throwExceptions(
false);
486 clp.addOutputSetupOptions(
true);
488 bool showTimers =
true;
490 clp.setOption(
"show-timers",
"no-show-timers", &showTimers,
"Determine if timers are shown or not" );
492 CommandLineProcessor::EParseCommandLineReturn
493 parse_return = clp.parse(argc,argv);
494 if( parse_return != CommandLineProcessor::PARSE_SUCCESSFUL )
498 out = VerboseObjectBase::getDefaultOStream();
502 result = masterTestComm<short int>(out);
503 if(!result) success =
false;
505 result = masterTestComm<int>(out);
506 if(!result) success =
false;
508 result = masterTestComm<long int>(out);
509 if(!result) success =
false;
514 ,out->getOutputToRootOnly() < 0
519 *out <<
"\nEnd Result: TEST PASSED\n";
524 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...