32#ifndef ALL_LB_HEADER_INCLUDED
33#define ALL_LB_HEADER_INCLUDED
44template <
class T,
class W>
class LB {
53 minSize = std::vector<T>(dim, 0.0);
61 virtual ~LB() =
default;
92 __FILE__, __func__, __LINE__,
93 "Currently only three dimensional vertices are supported.");
121 virtual void setWork(
const std::vector<W> &w) {
122 for(
const W& val : w) {
123 if(std::isnan(val)){
throw InvalidArgumentException(__FILE__, __func__, __LINE__,
"NAN is not a valid value for work.");}
124 if(std::isinf(val)){
throw InvalidArgumentException(__FILE__, __func__, __LINE__,
"inf is not a valid value for work.");}
199 double localSum = 0.0;
200 for (
auto i = this->
work.begin(); i < this->
work.end(); ++i)
206 MPI_Allreduce(&localSum, &globalMin, 1, MPI_DOUBLE, MPI_MIN, this->
globalComm);
207 MPI_Allreduce(&localSum, &globalMax, 1, MPI_DOUBLE, MPI_MAX, this->
globalComm);
208 return (1.0 - (
double)(globalMax - globalMin) /
209 (
double)(globalMax + globalMin));
virtual void setMinDomainSize(const std::vector< T > &minSize)
LB(const int dim, const T g)
virtual void setVertices(const std::vector< Point< T > > &vertices_in)
void resizeVertices(const int newSize)
virtual void setCommunicator(const MPI_Comm comm)
virtual void setDimension(const int d)
virtual void setWork(const W w)
virtual void setAdditionalData(const void *data)=0
int localRank
local rank in the used MPI communicator
std::vector< W > work
local work
MPI_Comm globalComm
used MPI communicator
int dimension
dimension of the used vertices
std::vector< T > sysSize
(orthogonal) system size
virtual std::vector< Point< T > > & getPrevVertices()
std::vector< T > & getNeighborVertices()
virtual void setGammaAutoCalc(const bool b)
std::vector< int > periodicity
periodicity of the MPI communicator / system
std::vector< int > global_dims
dimensions of the global process grid
virtual void setup()=0
abstract definition of the setup method
virtual const T getGamma()
virtual void setWork(const std::vector< W > &w)
virtual void setGamma(const T g)
virtual void setSysSize(const std::vector< T > &newSysSize)
virtual const bool getGammaAutoCalc()
bool gammaAutoCalc
turn on/off autocalculation of gamma factor if present
virtual W getEstimatedEfficiency()=0
virtual std::vector< W > & getWork()
virtual int getDimension()
std::vector< T > neighborVertices
vertices describing neighboring domains
virtual const std::vector< T > & getSysSize()
virtual ~LB()=default
destructor
std::vector< int > local_coords
cartesian coordinates of the local domain in the process grid
virtual std::vector< Point< T > > & getVertices()
std::vector< Point< T > > prevVertices
original vertices before previous balancing step
virtual const std::vector< T > & getMinDomainSize()
virtual void balance(const int step)=0
abstract definition of the balancing method
virtual std::vector< int > & getNeighbors()=0
std::vector< Point< T > > vertices
local vertices after previous balancing step