walberla::timing::TimingPool< TP > Class Template Reference

#include <TimingPool.h>

Public Member Functions

Construction & Destruction
 TimingPool ()=default
 
Reduction
shared_ptr< TimingPool< TP > > getReduced (ReduceType rt=REDUCE_TOTAL, int targetWorldRank=0) const
 Returns a reduced TimingPool, holding information from all processes. More...
 
Utility
void print (std::ostream &os) const
 Prints the timing pool to the given ostream. More...
 
void printMatlab (std::ostream &os) const
 Prints the TimingPool in a format that can be read by MATLAB for plotting. More...
 
void printCSV (std::ostream &os, const std::string &separator=", ") const
 Prints the TimingPool in a format that can be read by MATLAB for plotting. More...
 
void unifyRegisteredTimersAcrossProcesses ()
 Ensures that on each process the same timers (= timers with identical names) are registered. More...
 
void logResultOnRoot (ReduceType rt=REDUCE_TOTAL, const bool unifyRegisteredTimers=false)
 Reduces the timing pool, and logs reduced result on root. More...
 
void merge (const TimingPool< TP > &tpToMerge, bool mergeDuplicates=true)
 Merges a second TimingPool into this TimingPool. More...
 
void clear ()
 Resets all timers. More...
 
bool empty () const
 

Friends

template<typename T , typename G , typename TP2 >
mpi::GenericSendBuffer< T, G > & operator<< (mpi::GenericSendBuffer< T, G > &buf, const TimingPool< TP2 > &tp)
 
template<typename T , typename TP2 >
mpi::GenericRecvBuffer< T > & operator>> (mpi::GenericRecvBuffer< T > &buf, TimingPool< TP2 > &tp)
 

Access to Timer Objects

using iterator = typename std::map< std::string, Timer< TP > >::iterator
 
using const_iterator = typename std::map< std::string, Timer< TP > >::const_iterator
 
iterator begin ()
 
const_iterator begin () const
 
iterator end ()
 
const_iterator end () const
 
Timer< TP > & operator[] (const std::string &n)
 
const Timer< TP > & operator[] (const std::string &n) const
 
bool timerExists (const std::string &n) const
 
bool registerTimer (const std::string &n)
 
ScopeTimer< TP > getScopeTimer (const std::string &n)
 

Reduction Helper

std::map< std::string, Timer< TP > > timerMap_
 
void mpiReduce (std::vector< double > &min, std::vector< double > &max, std::vector< double > &val, std::vector< double > &valSq, std::vector< uint32_t > &count, int targetRank, shared_ptr< TimingPool< TP > > &out) const
 Executes an mpi reduction of the given vectors, and stores the result in a new TimingPool. More...
 

Member Typedef Documentation

◆ const_iterator

template<typename TP >
using walberla::timing::TimingPool< TP >::const_iterator = typename std::map<std::string, Timer<TP> >::const_iterator

◆ iterator

template<typename TP >
using walberla::timing::TimingPool< TP >::iterator = typename std::map<std::string, Timer<TP> >::iterator

Constructor & Destructor Documentation

◆ TimingPool()

template<typename TP >
walberla::timing::TimingPool< TP >::TimingPool ( )
default

Member Function Documentation

◆ begin() [1/2]

template<typename TP >
iterator walberla::timing::TimingPool< TP >::begin ( )
inline

◆ begin() [2/2]

template<typename TP >
const_iterator walberla::timing::TimingPool< TP >::begin ( ) const
inline

◆ clear()

template<typename TP >
void walberla::timing::TimingPool< TP >::clear

Resets all timers.

◆ empty()

template<typename TP >
bool walberla::timing::TimingPool< TP >::empty ( ) const
inline

◆ end() [1/2]

template<typename TP >
iterator walberla::timing::TimingPool< TP >::end ( )
inline

◆ end() [2/2]

template<typename TP >
const_iterator walberla::timing::TimingPool< TP >::end ( ) const
inline

◆ getReduced()

template<typename TP >
shared_ptr< TimingPool< TP > > walberla::timing::TimingPool< TP >::getReduced ( ReduceType  rt = REDUCE_TOTAL,
int  targetRank = 0 
) const

Returns a reduced TimingPool, holding information from all processes.

This function has to be called on all processes. The TimingPool has to have the same number of timers on each processor, and the timers have to have the same name.

Parameters
rtSpecified the method how the reduction is done. See documentation for ReduceType
targetRankthe world rank of the target process. Or negative value for an all-reduction operation
Returns
a nonzero shared pointer to a timing pool on the process with rank "targetRank" and a zero pointer otherwise. If targetRank < 0, a valid pointer is returned on all processes (MPI_Allreduce operation)

◆ getScopeTimer()

template<typename TP >
ScopeTimer< TP > walberla::timing::TimingPool< TP >::getScopeTimer ( const std::string &  n)
inline

◆ logResultOnRoot()

template<typename TP >
void walberla::timing::TimingPool< TP >::logResultOnRoot ( ReduceType  rt = REDUCE_TOTAL,
const bool  unifyRegisteredTimers = false 
)

Reduces the timing pool, and logs reduced result on root.

Has to be called by all processes, not just on root process.

◆ merge()

template<typename TP >
void walberla::timing::TimingPool< TP >::merge ( const TimingPool< TP > &  tpToMerge,
bool  mergeDuplicates = true 
)

Merges a second TimingPool into this TimingPool.

Parameters
tpToMergethe second TimingPool, which is merged
mergeDuplicatesparameter has an effect when both TimingPool's have Timers with same name if true: the timers with the same name are merged together if false: the timer in this TimingPool is not changed, the content of the second timer with the same name is neglected

◆ mpiReduce()

template<typename TP >
void walberla::timing::TimingPool< TP >::mpiReduce ( std::vector< double > &  min,
std::vector< double > &  max,
std::vector< double > &  val,
std::vector< double > &  valSq,
std::vector< uint32_t > &  count,
int  targetRank,
shared_ptr< TimingPool< TP > > &  out 
) const
protected

Executes an mpi reduction of the given vectors, and stores the result in a new TimingPool.

◆ operator[]() [1/2]

template<typename TP >
Timer< TP > & walberla::timing::TimingPool< TP >::operator[] ( const std::string &  n)
inline

◆ operator[]() [2/2]

template<typename TP >
const Timer< TP > & walberla::timing::TimingPool< TP >::operator[] ( const std::string &  n) const
inline

◆ print()

template<typename TP >
void walberla::timing::TimingPool< TP >::print ( std::ostream &  os) const

Prints the timing pool to the given ostream.

No reduction is done!

If the result of all processes should be printed, first a reduction has to be done:

tp.print(cout); // prints only the local time measurements
shared_ptr< WcTimingPool> reduced = tp.getReduced(REDUCE_TOTAL, 0);
reduced->print(cout);
}

◆ printCSV()

template<typename TP >
void walberla::timing::TimingPool< TP >::printCSV ( std::ostream &  os,
const std::string &  separator = ", " 
) const

Prints the TimingPool in a format that can be read by MATLAB for plotting.

◆ printMatlab()

template<typename TP >
void walberla::timing::TimingPool< TP >::printMatlab ( std::ostream &  os) const

Prints the TimingPool in a format that can be read by MATLAB for plotting.

◆ registerTimer()

template<typename TP >
bool walberla::timing::TimingPool< TP >::registerTimer ( const std::string &  n)
inline

◆ timerExists()

template<typename TP >
bool walberla::timing::TimingPool< TP >::timerExists ( const std::string &  n) const
inline

◆ unifyRegisteredTimersAcrossProcesses()

template<typename TP >
void walberla::timing::TimingPool< TP >::unifyRegisteredTimersAcrossProcesses

Ensures that on each process the same timers (= timers with identical names) are registered.

Has to be called by all processes, not just on root process.

Friends And Related Function Documentation

◆ operator<<

template<typename TP >
template<typename T , typename G , typename TP2 >
mpi::GenericSendBuffer<T,G>& operator<< ( mpi::GenericSendBuffer< T, G > &  buf,
const TimingPool< TP2 > &  tp 
)
friend

◆ operator>>

template<typename TP >
template<typename T , typename TP2 >
mpi::GenericRecvBuffer<T>& operator>> ( mpi::GenericRecvBuffer< T > &  buf,
TimingPool< TP2 > &  tp 
)
friend

Member Data Documentation

◆ timerMap_

template<typename TP >
std::map<std::string, Timer<TP> > walberla::timing::TimingPool< TP >::timerMap_
protected

The documentation for this class was generated from the following files:
#define WALBERLA_ROOT_SECTION()
Definition: MPIManager.h:174
@ REDUCE_TOTAL
Collects all timing samples from all processes and accumulates the data.
Definition: ReduceType.h:37