ALL 0.9.4
A Loadbalacing Library
Loading...
Searching...
No Matches
ALL_LB.hpp
Go to the documentation of this file.
1/*
2Copyright 2018-2026 Rene Halver, Forschungszentrum Juelich GmbH, Germany
3Copyright 2018-2026 Godehard Sutmann, Forschungszentrum Juelich GmbH, Germany
4Copyright 2026-2026 Jonas Kroschewski, Forschungszentrum Juelich GmbH, Germany
5
6Redistribution and use in source and binary forms, with or without modification,
7are permitted provided that the following conditions are met:
8
91. Redistributions of source code must retain the above copyright notice, this
10 list of conditions and the following disclaimer.
11
122. Redistributions in binary form must reproduce the above copyright notice,
13this list of conditions and the following disclaimer in the documentation and/or
14 other materials provided with the distribution.
15
163. Neither the name of the copyright holder nor the names of its contributors
17may be used to endorse or promote products derived from this software without
18 specific prior written permission.
19
20THIS SOFTWARE IS PROVIDED BY THE COPYRIGHT HOLDERS AND CONTRIBUTORS "AS IS" AND
21ANY EXPRESS OR IMPLIED WARRANTIES, INCLUDING, BUT NOT LIMITED TO, THE IMPLIED
22WARRANTIES OF MERCHANTABILITY AND FITNESS FOR A PARTICULAR PURPOSE ARE
23DISCLAIMED. IN NO EVENT SHALL THE COPYRIGHT HOLDER OR CONTRIBUTORS BE LIABLE FOR
24ANY DIRECT, INDIRECT, INCIDENTAL, SPECIAL, EXEMPLARY, OR CONSEQUENTIAL DAMAGES
25(INCLUDING, BUT NOT LIMITED TO, PROCUREMENT OF SUBSTITUTE GOODS OR SERVICES;
26LOSS OF USE, DATA, OR PROFITS; OR BUSINESS INTERRUPTION) HOWEVER CAUSED AND ON
27ANY THEORY OF LIABILITY, WHETHER IN CONTRACT, STRICT LIABILITY, OR TORT
28(INCLUDING NEGLIGENCE OR OTHERWISE) ARISING IN ANY WAY OUT OF THE USE OF THIS
29SOFTWARE, EVEN IF ADVISED OF THE POSSIBILITY OF SUCH DAMAGE.
30*/
31
32#ifndef ALL_LB_HEADER_INCLUDED
33#define ALL_LB_HEADER_INCLUDED
34
35#include "ALL_Point.hpp"
36#include <mpi.h>
37#include <vector>
38#include <iostream>
39
40namespace ALL {
41
44template <class T, class W> class LB {
45public:
50 LB(const int dim, const T g) : gamma(g) {
51 // set allowed minimum size to zero
52 setDimension(dim);
53 minSize = std::vector<T>(dim, 0.0);
54 periodicity.resize(dim);
55 local_coords.resize(dim);
56 global_dims.resize(dim);
57 neighborVertices.resize(0);
58 }
59
61 virtual ~LB() = default;
62
64 virtual void setup() = 0;
65
67 virtual void balance(const int step) = 0;
68
72 virtual std::vector<int> &getNeighbors() = 0;
73
77 virtual void setVertices(const std::vector<Point<T>> &vertices_in) {
78 // as a basic assumption the number of resulting vertices
79 // is equal to the number of input vertices:
80 // exceptions:
81 // VORONOI
82 vertices = vertices_in;
83 prevVertices = vertices_in;
84 }
85
89 virtual void setDimension(const int d) {
90 if (d != 3)
92 __FILE__, __func__, __LINE__,
93 "Currently only three dimensional vertices are supported.");
94 else
95 dimension = d;
96 }
97
100 virtual int getDimension() { return dimension; }
101
104 virtual void setGamma(const T g) { gamma = g; }
105
108 virtual const T getGamma() { return gamma; }
109
110 // method to set the automatic calculation of the correction value
111 // @param[in] b the boolean to activate/deactivate automatic calculation
112 virtual void setGammaAutoCalc(const bool b) { gammaAutoCalc = b; }
113
116 virtual const bool getGammaAutoCalc() { return gammaAutoCalc; }
117
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.");}
125 }
126 work = w;
127 }
128
131 virtual void setWork(const W w) {
132 work.resize(1);
133 work.at(0) = w;
134 }
135
138 virtual void setMinDomainSize(const std::vector<T> &minSize) {
139 this->minSize = minSize;
140 }
141
144 virtual const std::vector<T> &getMinDomainSize() { return minSize; }
145
148 virtual void setSysSize(const std::vector<T> &newSysSize) {
149 sysSize = newSysSize;
150 }
151
154 virtual const std::vector<T> &getSysSize() { return sysSize; }
155
159 virtual std::vector<W> &getWork() { return work; }
160
163 virtual std::vector<Point<T>> &getVertices() { return vertices; }
164
170 virtual std::vector<Point<T>> &getPrevVertices() { return prevVertices; };
171
174 virtual void setCommunicator(const MPI_Comm comm) {
175 // store communicator
176 globalComm = comm;
177 // get local rank within communicator
178 MPI_Comm_rank(comm, &localRank);
179 }
180
183 int getNVertices() { return vertices.size(); }
184
188 virtual void setAdditionalData(const void *data) = 0;
189
193 virtual W getEstimatedEfficiency() = 0;
194
198 {
199 double localSum = 0.0;
200 for (auto i = this->work.begin(); i < this->work.end(); ++i)
201 {
202 localSum += *i;
203 }
204 double globalMin;
205 double globalMax;
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));
210 }
211
212
216 std::vector<T> &getNeighborVertices() {
217 return neighborVertices; }
218
219protected:
227 MPI_Comm globalComm;
232 std::vector<T> minSize;
234 std::vector<T> sysSize;
236 std::vector<W> work;
238 std::vector<Point<T>> vertices;
240 std::vector<Point<T>> prevVertices;
242 std::vector<int> global_dims;
244 std::vector<int> local_coords;
246 std::vector<int> periodicity;
248 std::vector<T> neighborVertices;
251 void resizeVertices(const int newSize) { vertices.resize(newSize); }
252
253};
254
255}//namespace ALL
256
257#endif // ALL_LB_HEADER_INCLUDED
virtual void setMinDomainSize(const std::vector< T > &minSize)
Definition ALL_LB.hpp:138
LB(const int dim, const T g)
Definition ALL_LB.hpp:50
virtual void setVertices(const std::vector< Point< T > > &vertices_in)
Definition ALL_LB.hpp:77
void resizeVertices(const int newSize)
Definition ALL_LB.hpp:251
virtual void setCommunicator(const MPI_Comm comm)
Definition ALL_LB.hpp:174
virtual void setDimension(const int d)
Definition ALL_LB.hpp:89
virtual void setWork(const W w)
Definition ALL_LB.hpp:131
virtual void setAdditionalData(const void *data)=0
int localRank
local rank in the used MPI communicator
Definition ALL_LB.hpp:229
std::vector< W > work
local work
Definition ALL_LB.hpp:236
std::vector< T > minSize
Definition ALL_LB.hpp:232
MPI_Comm globalComm
used MPI communicator
Definition ALL_LB.hpp:227
int dimension
dimension of the used vertices
Definition ALL_LB.hpp:225
std::vector< T > sysSize
(orthogonal) system size
Definition ALL_LB.hpp:234
virtual std::vector< Point< T > > & getPrevVertices()
Definition ALL_LB.hpp:170
int getNVertices()
Definition ALL_LB.hpp:183
std::vector< T > & getNeighborVertices()
Definition ALL_LB.hpp:216
virtual void setGammaAutoCalc(const bool b)
Definition ALL_LB.hpp:112
std::vector< int > periodicity
periodicity of the MPI communicator / system
Definition ALL_LB.hpp:246
std::vector< int > global_dims
dimensions of the global process grid
Definition ALL_LB.hpp:242
virtual void setup()=0
abstract definition of the setup method
virtual const T getGamma()
Definition ALL_LB.hpp:108
virtual void setWork(const std::vector< W > &w)
Definition ALL_LB.hpp:121
virtual void setGamma(const T g)
Definition ALL_LB.hpp:104
virtual void setSysSize(const std::vector< T > &newSysSize)
Definition ALL_LB.hpp:148
double getEfficiency()
Definition ALL_LB.hpp:197
virtual const bool getGammaAutoCalc()
Definition ALL_LB.hpp:116
bool gammaAutoCalc
turn on/off autocalculation of gamma factor if present
Definition ALL_LB.hpp:223
virtual W getEstimatedEfficiency()=0
virtual std::vector< W > & getWork()
Definition ALL_LB.hpp:159
T gamma
correction factor
Definition ALL_LB.hpp:221
virtual int getDimension()
Definition ALL_LB.hpp:100
std::vector< T > neighborVertices
vertices describing neighboring domains
Definition ALL_LB.hpp:248
virtual const std::vector< T > & getSysSize()
Definition ALL_LB.hpp:154
virtual ~LB()=default
destructor
std::vector< int > local_coords
cartesian coordinates of the local domain in the process grid
Definition ALL_LB.hpp:244
virtual std::vector< Point< T > > & getVertices()
Definition ALL_LB.hpp:163
std::vector< Point< T > > prevVertices
original vertices before previous balancing step
Definition ALL_LB.hpp:240
virtual const std::vector< T > & getMinDomainSize()
Definition ALL_LB.hpp:144
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
Definition ALL_LB.hpp:238
Definition ALL.hpp:78