Skip to content
Snippets Groups Projects
Commit 16b0f63c authored by Sebastian Eibl's avatar Sebastian Eibl
Browse files

added reduction to Timer

parent 6d720b5c
Branches
No related merge requests found
......@@ -25,10 +25,12 @@
#pragma once
#include "CpuPolicy.h"
#include "ReduceType.h"
#include "WcPolicy.h"
#include "core/DataTypes.h"
#include "core/mpi/RecvBuffer.h"
#include "core/mpi/Reduce.h"
#include "core/mpi/SendBuffer.h"
#include <iomanip>
......@@ -155,7 +157,7 @@ public:
//@}
//*******************************************************************************************************************
private:
private:
uint_t counter_; //!< Number of performed time measurements.
double start_; //!< Start of the current time measurement.
......@@ -211,14 +213,14 @@ inline Timer<TP>::Timer()
template< typename TP>
inline Timer<TP>::Timer( uint_t _counter, double _min, double _max,
double _total, double _sumOfSquares )
: counter_ ( _counter )
, start_ ( 0.0 )
, end_ ( 0.0 )
, time_ ( _total )
, sumOfSquares_( _sumOfSquares )
, min_ ( _min )
, max_ ( _max )
, last_ ( 0.0 )
: counter_ ( _counter )
, start_ ( 0.0 )
, end_ ( 0.0 )
, time_ ( _total )
, sumOfSquares_( _sumOfSquares )
, min_ ( _min )
, max_ ( _max )
, last_ ( 0.0 )
{
}
//**********************************************************************************************************************
......@@ -455,6 +457,80 @@ inline void Timer<TP>::merge( const Timer<TP> & other )
//**********************************************************************************************************************
//======================================================================================================================
//
// REDUCTION
//
//======================================================================================================================
//**********************************************************************************************************************
/*! Returns a reduced Timer, holding information from all processes
*
* \param timer Timer which should be reduced
* \param rt Specified the method how the reduction is done. See documentation for ReduceType
* \param targetRank the world rank of the target process. Or negative value for an all-reduction
* operation
*
* \return a nonzero shared pointer to a Timer on the process with rank "targetRank"
* and a zero pointer otherwise. For the all-reduction a valid timer is returned on all processes.
**********************************************************************************************************************/
template< typename TP >
shared_ptr<Timer<TP> > getReduced( Timer<TP>& timer, ReduceType rt, int targetRank )
{
WALBERLA_NON_MPI_SECTION() {
return make_shared<Timer<TP> >( timer );
}
double val; //value to be reduced
switch (rt)
{
case REDUCE_MIN :
val = timer.min();
break;
case REDUCE_AVG :
val = timer.average();
break;
case REDUCE_MAX :
val = timer.max();
break;
case REDUCE_TOTAL:
val = timer.total();
break;
default:
WALBERLA_ABORT( "Unknown reduce type" );
break;
}
double min;
double max;
double total;
double sumOfSquares;
if (targetRank >= 0)
{
min = mpi::reduce(val, mpi::MIN, targetRank);
max = mpi::reduce(val, mpi::MAX, targetRank);
total = mpi::reduce(val, mpi::SUM, targetRank);
sumOfSquares = mpi::reduce(val*val, mpi::SUM, targetRank);
} else
{
min = mpi::allReduce(val, mpi::MIN);
max = mpi::allReduce(val, mpi::MAX);
total = mpi::allReduce(val, mpi::SUM);
sumOfSquares = mpi::allReduce(val*val, mpi::SUM);
}
//uint_t counter, double min, double max, double total, double sumOfSquares
if ( targetRank < 0 || targetRank == MPIManager::instance()->worldRank() )
return make_shared<Timer<TP> >( mpi::MPIManager::instance()->numProcesses(), min, max, total, sumOfSquares );
return nullptr;
}
//======================================================================================================================
//
// OSTREAM OVERLOAD
......@@ -478,38 +554,38 @@ std::ostream & operator<< ( std::ostream & os, const Timer<TP> & timer )
//
//======================================================================================================================
template< typename T, // Element type of SendBuffer
typename G, // Growth policy of SendBuffer
typename TP > // Element type of vector
mpi::GenericSendBuffer<T,G>& operator<<( mpi::GenericSendBuffer<T,G> & buf, const Timer<TP> & t )
{
buf.addDebugMarker( "ti" );
buf << t.counter_;
buf << t.start_;
buf << t.end_;
buf << t.time_;
buf << t.sumOfSquares_;
buf << t.min_;
buf << t.max_;
buf << t.last_;
return buf;
}
template< typename T, // Element type of SendBuffer
typename G, // Growth policy of SendBuffer
typename TP > // Element type of vector
mpi::GenericSendBuffer<T,G>& operator<<( mpi::GenericSendBuffer<T,G> & buf, const Timer<TP> & t )
{
buf.addDebugMarker( "ti" );
buf << t.counter_;
buf << t.start_;
buf << t.end_;
buf << t.time_;
buf << t.sumOfSquares_;
buf << t.min_;
buf << t.max_;
buf << t.last_;
return buf;
}
template< typename T, // Element type of RecvBuffer
typename TP > // Element type of vector
mpi::GenericRecvBuffer<T>& operator>>( mpi::GenericRecvBuffer<T> & buf, Timer<TP> & t )
{
buf.readDebugMarker( "ti" );
buf >> t.counter_;
buf >> t.start_;
buf >> t.end_;
buf >> t.time_;
buf >> t.sumOfSquares_;
buf >> t.min_;
buf >> t.max_;
buf >> t.last_;
return buf;
}
template< typename T, // Element type of RecvBuffer
typename TP > // Element type of vector
mpi::GenericRecvBuffer<T>& operator>>( mpi::GenericRecvBuffer<T> & buf, Timer<TP> & t )
{
buf.readDebugMarker( "ti" );
buf >> t.counter_;
buf >> t.start_;
buf >> t.end_;
buf >> t.time_;
buf >> t.sumOfSquares_;
buf >> t.min_;
buf >> t.max_;
buf >> t.last_;
return buf;
}
} //namespace timing
......
......@@ -159,6 +159,9 @@ waLBerla_execute_test( NAME SetSelectableObjectTest )
# timing #
##########
waLBerla_compile_test( FILES timing/ParallelTimerTest.cpp )
waLBerla_execute_test( NAME ParallelTimerTest PROCESSES 2 )
waLBerla_compile_test( FILES timing/TimerTest.cpp )
waLBerla_execute_test( NAME TimerTest )
......
//======================================================================================================================
//
// This file is part of waLBerla. waLBerla is free software: you can
// redistribute it and/or modify it under the terms of the GNU General Public
// License as published by the Free Software Foundation, either version 3 of
// the License, or (at your option) any later version.
//
// waLBerla is distributed in the hope that it will be useful, but WITHOUT
// ANY WARRANTY; without even the implied warranty of MERCHANTABILITY or
// FITNESS FOR A PARTICULAR PURPOSE. See the GNU General Public License
// for more details.
//
// You should have received a copy of the GNU General Public License along
// with waLBerla (see COPYING.txt). If not, see <http://www.gnu.org/licenses/>.
//
//! \file ParallelTimerTest.cpp
//! \author Sebastian Eibl <sebastian.eibl@web.de>
//
//======================================================================================================================
#include "core/debug/TestSubsystem.h"
#include "core/Environment.h"
#include "core/math/Constants.h"
#include "core/timing/StaticPolicy.h"
#include "core/timing/Timer.h"
#include <cmath>
namespace walberla {
void reductionTest()
{
using namespace timing;
StaticPolicy::setTime(0.0);
Timer<StaticPolicy> t0;
t0.start();
StaticPolicy::setTime(2.0);
t0.end();
t0.start();
StaticPolicy::setTime(6.0);
t0.end();
WALBERLA_CHECK_FLOAT_EQUAL( t0.min(), 2.0 );
WALBERLA_CHECK_FLOAT_EQUAL( t0.max(), 4.0 );
WALBERLA_CHECK_FLOAT_EQUAL( t0.average(), 3.0 );
WALBERLA_CHECK_EQUAL( t0.getCounter(), 2 );
auto timer_reduced = getReduced(t0, REDUCE_TOTAL, -1); //allreduce
WALBERLA_CHECK_FLOAT_EQUAL( timer_reduced->min(), 6.0 );
WALBERLA_CHECK_FLOAT_EQUAL( timer_reduced->max(), 6.0 );
WALBERLA_CHECK_FLOAT_EQUAL( timer_reduced->average(), 6.0 );
WALBERLA_CHECK_EQUAL( timer_reduced->getCounter(), 2 );
timer_reduced = getReduced(t0, REDUCE_TOTAL, 0);
if (mpi::MPIManager::instance()->worldRank() == 0)
{
WALBERLA_CHECK_FLOAT_EQUAL( timer_reduced->min(), 6.0 );
WALBERLA_CHECK_FLOAT_EQUAL( timer_reduced->max(), 6.0 );
WALBERLA_CHECK_FLOAT_EQUAL( timer_reduced->average(), 6.0 );
WALBERLA_CHECK_EQUAL( timer_reduced->getCounter(), 2 );
} else
{
WALBERLA_CHECK_EQUAL( timer_reduced, nullptr );
}
}
}
int main( int argc, char ** argv )
{
walberla::debug::enterTestMode();
walberla::Environment env(argc, argv);
WALBERLA_UNUSED(env);
walberla::reductionTest();
}
0% or .
You are about to add 0 people to the discussion. Proceed with caution.
Finish editing this message first!
Please register or to comment