Skip to content
Snippets Groups Projects
Commit db2c3599 authored by Tomáš Oberhuber's avatar Tomáš Oberhuber Committed by Tomáš Oberhuber
Browse files

String send and receive was changed from methods to functions.

parent 0c1362a9
No related branches found
No related tags found
1 merge request!29Revision
......@@ -29,15 +29,23 @@ ADD_EXECUTABLE( MathExample MathExample.cpp )
ADD_EXECUTABLE( ParameterContainerExample ParameterContainerExample.cpp )
ADD_EXECUTABLE( StringExample StringExample.cpp )
ADD_CUSTOM_COMMAND( OUTPUT StringExample.out
COMMAND StringExample
DEPENDS StringExample
COMMENT Executing StringExample )
EXECUTE_PROCESS( COMMAND StringExample OUTPUT_FILE StringExample.out )
ADD_EXECUTABLE( StringExampleGetAllocatedSize StringExampleGetAllocatedSize.cpp )
EXECUTE_PROCESS( COMMAND StringExampleGetAllocatedSize OUTPUT_FILE StringExampleGetAllocatedSize.out )
ADD_EXECUTABLE( StringExampleReplace StringExampleReplace.cpp )
EXECUTE_PROCESS( COMMAND StringExampleReplace OUTPUT_FILE StringExampleReplace.out )
ADD_EXECUTABLE( StringExampleSetSize StringExampleSetSize.cpp )
EXECUTE_PROCESS( COMMAND StringExampleSetSize OUTPUT_FILE StringExampleSetSize.out )
ADD_EXECUTABLE( StringExampleSplit StringExampleSplit.cpp )
EXECUTE_PROCESS( COMMAND StringExampleSpli OUTPUT_FILE StringExampleSplit.out )
ADD_EXECUTABLE( StringExampleStrip StringExampleStrip.cpp )
EXECUTE_PROCESS( COMMAND StringExampleStrip OUTPUT_FILE StringExampleStrip.out )
ADD_EXECUTABLE( TimerExample TimerExample.cpp )
ADD_EXECUTABLE( VectorExample VectorExample.cpp )
......@@ -25,7 +25,7 @@ else
__tnl_mpi_print_stream_ << "Node " << TNL::Communicators::MpiCommunicator::GetRank() << " of " \
<< TNL::Communicators::MpiCommunicator::GetSize() << " : " << message << std::endl; \
TNL::String __tnl_mpi_print_string_( __tnl_mpi_print_stream_.str().c_str() ); \
__tnl_mpi_print_string_.send( 0, std::numeric_limits< int >::max() ); \
send( __tnl_mpi_print_string_, 0, std::numeric_limits< int >::max() ); \
} \
else \
{ \
......@@ -35,7 +35,7 @@ else
__tnl_mpi_print_j++ ) \
{ \
TNL::String __tnl_mpi_print_string_; \
__tnl_mpi_print_string_.receive( __tnl_mpi_print_j, std::numeric_limits< int >::max() ); \
receive( __tnl_mpi_print_string_, __tnl_mpi_print_j, std::numeric_limits< int >::max() ); \
std::cerr << __tnl_mpi_print_string_; \
} \
} \
......@@ -78,7 +78,7 @@ else
__tnl_mpi_print_stream_ << "Node " << TNL::Communicators::MpiCommunicator::GetRank() << " of " \
<< TNL::Communicators::MpiCommunicator::GetSize() << " : " << message << std::endl; \
TNL::String __tnl_mpi_print_string_( __tnl_mpi_print_stream_.str().c_str() ); \
__tnl_mpi_print_string_.send( 0, std::numeric_limits< int >::max() ); \
send( __tnl_mpi_print_string_, 0, std::numeric_limits< int >::max() ); \
} \
} \
else \
......@@ -94,7 +94,7 @@ else
if( __tnl_mpi_print_cond ) \
{ \
TNL::String __tnl_mpi_print_string_; \
__tnl_mpi_print_string_.receive( __tnl_mpi_print_j, std::numeric_limits< int >::max() ); \
receive( __tnl_mpi_print_string_, __tnl_mpi_print_j, std::numeric_limits< int >::max() ); \
std::cerr << __tnl_mpi_print_string_; \
} \
} \
......
......@@ -337,36 +337,33 @@ class String
* \include StringExampleSplit.out
*/
std::vector< String > split( const char separator = ' ', SplitSkipEmpty skipEmpty = DoNotSkipEmpty ) const;
#ifdef HAVE_MPI
/**
* \brief Sends the string to the target MPI process.
*/
void send( int target, int tag = 0, MPI_Comm mpi_comm = MPI_COMM_WORLD );
/**
* \brief Receives a string from the source MPI process.
*/
void receive( int source, int tag = 0, MPI_Comm mpi_comm = MPI_COMM_WORLD );
//! Broadcast to other nodes in MPI cluster
// void MPIBcast( int root, MPI_Comm mpi_comm = MPI_COMM_WORLD );
#endif
};
/// \brief Returns concatenation of \e string1 and \e string2.
/**
* \brief Returns concatenation of \e string1 and \e string2.
*/
String operator+( char string1, const String& string2 );
/// \brief Returns concatenation of \e string1 and \e string2.
/**
* \brief Returns concatenation of \e string1 and \e string2.
*/
String operator+( const char* string1, const String& string2 );
/// \brief Returns concatenation of \e string1 and \e string2.
/**
* \brief Returns concatenation of \e string1 and \e string2.
*/
String operator+( const std::string& string1, const String& string2 );
/// \brief Performs the string output to a stream
/**
* \brief Writes the string \e str to given \e stream
*/
std::ostream& operator<<( std::ostream& stream, const String& str );
/**
* \brief Converts \e value of type \e T to a String.
*
* \tparam T can be any type fir which operator << is defined.
*/
template< typename T >
String convertToString( const T& value )
{
......@@ -375,12 +372,48 @@ String convertToString( const T& value )
return String( str.str().data() );
}
/**
* \brief Specialization of function \ref conertToString for boolean.
*
* The boolean type is converted to 'true' ot 'false'.
*/
template<> inline String convertToString( const bool& b )
{
if( b ) return "true";
return "false";
}
#ifdef HAVE_MPI
/**
* \brief Sends the string to the target MPI process.
*
* @param str string to be sent
* @param target target MPI process ID
* @param tag MPI tag
* @param mpi_comm MPI communication group
*/
void send( const String& str, int target, int tag = 0, MPI_Comm mpi_comm = MPI_COMM_WORLD );
/**
* \brief Receives a string from the source MPI process.
*/
/**
* \brief Receives a string from the target MPI process.
*
* @param str says where the received string is to be saved to
* @param source source MPI process ID
* @param tag MPI tag
* @param mpi_comm MPI communication group
*/
void receive( String& str, int source, int tag = 0, MPI_Comm mpi_comm = MPI_COMM_WORLD );
//! Broadcast to other nodes in MPI cluster
// void MPIBcast( String& str, int root, MPI_Comm mpi_comm = MPI_COMM_WORLD );
#endif
} // namespace TNL
#include <TNL/String.hpp>
......@@ -228,27 +228,47 @@ String::split( const char separator, SplitSkipEmpty skipEmpty ) const
return parts;
}
inline String operator+( char string1, const String& string2 )
{
return convertToString( string1 ) + string2;
}
inline String operator+( const char* string1, const String& string2 )
{
return String( string1 ) + string2;
}
inline String operator+( const std::string& string1, const String& string2 )
{
return String( string1 ) + string2;
}
inline std::ostream& operator<<( std::ostream& stream, const String& str )
{
stream << str.getString();
return stream;
}
#ifdef HAVE_MPI
inline void String::send( int target, int tag, MPI_Comm mpi_comm )
inline void send( const String& str, int target, int tag, MPI_Comm mpi_comm )
{
int size = this->getSize();
int size = str.getSize();
MPI_Send( &size, 1, MPI_INT, target, tag, mpi_comm );
MPI_Send( this->getString(), this->length(), MPI_CHAR, target, tag, mpi_comm );
MPI_Send( str.getString(), str.length(), MPI_CHAR, target, tag, mpi_comm );
}
inline void String::receive( int source, int tag, MPI_Comm mpi_comm )
inline void receive( String& str, int source, int tag, MPI_Comm mpi_comm )
{
int size;
MPI_Status status;
MPI_Recv( &size, 1, MPI_INT, source, tag, mpi_comm, &status );
this->setSize( size );
MPI_Recv( const_cast< void* >( ( const void* ) this->data() ), size, MPI_CHAR, source, tag, mpi_comm, &status );
str.setSize( size );
MPI_Recv( const_cast< void* >( ( const void* ) str.data() ), size, MPI_CHAR, source, tag, mpi_comm, &status );
}
/*
inline void String :: MPIBcast( int root, MPI_Comm comm )
{
#ifdef USE_MPI
int iproc;
MPI_Comm_rank( MPI_COMM_WORLD, &iproc );
TNL_ASSERT( string, );
......@@ -265,30 +285,9 @@ inline void String :: MPIBcast( int root, MPI_Comm comm )
}
MPI_Bcast( string, len + 1, MPI_CHAR, root, comm );
#endif
}
*/
#endif
inline String operator+( char string1, const String& string2 )
{
return convertToString( string1 ) + string2;
}
inline String operator+( const char* string1, const String& string2 )
{
return String( string1 ) + string2;
}
inline String operator+( const std::string& string1, const String& string2 )
{
return String( string1 ) + string2;
}
inline std::ostream& operator<<( std::ostream& stream, const String& str )
{
stream << str.getString();
return stream;
}
} // namespace TNL
0% Loading or .
You are about to add 0 people to the discussion. Proceed with caution.
Please register or to comment