SetLoad.cpp
1 #include "SetLoad.hpp"
2 #ifdef RICH_MPI
3 #include <mpi.h>
4 #endif
5 
6 
7 #ifdef RICH_MPI
8 
9 namespace
10 {
11  double GetLoad(Tessellation const& tess)
12  {
13  int ws;
14  MPI_Comm_size(MPI_COMM_WORLD, &ws);
15  int n=tess.GetPointNo();
16  int result=10;
17  int total=n;
18  MPI_Reduce(&n, &result, 1, MPI_INT, MPI_MAX, 0, MPI_COMM_WORLD);
19  MPI_Reduce(&n, &total, 1, MPI_INT, MPI_SUM, 0, MPI_COMM_WORLD);
20  return result*ws*1.0/(total*1.0);
21  }
22 
23  void RunLoadBalance(Tessellation &tproc,vector<Vector2D> &points,OuterBoundary const&
24  outer,int Niter,double tload,double speed,int mode,bool Rmin)
25  {
26  int ws;
27  MPI_Comm_size(MPI_COMM_WORLD, &ws);
28  int rank = 0;
29  MPI_Comm_rank(MPI_COMM_WORLD, &rank);
30  double round;
31  if(mode==1)
32  round=2;
33  else
34  round=1.6;
35  vector<size_t> selfindex;
36  vector<vector<int> > sentpoints;
37  vector<int> sentproc;
38  ConstNumberPerProc procmove(outer,speed,round,mode,Rmin);
39  VoronoiMesh local;
40  local.GetMeshPoints() = points;
41  local.SetPointNo(static_cast<int>(points.size()));
42  for(int i=0;i<Niter;++i)
43  {
44  MPI_Barrier(MPI_COMM_WORLD);
45  procmove.Update(tproc,local);
46  points = local.UpdateMPIPoints(tproc, rank, points, &outer, selfindex, sentproc, sentpoints);
47  local.GetMeshPoints() = points;
48  local.SetPointNo(static_cast<int>(points.size()));
49  MPI_Barrier(MPI_COMM_WORLD);
50  double load=GetLoad(local);
51  MPI_Bcast(&load, 1, MPI_DOUBLE, 0, MPI_COMM_WORLD);
52  if(load<tload)
53  {
54  break;
55  }
56  }
57  }
58 }
59 
60 void SetLoad(Tessellation &tproc,vector<Vector2D> &points,OuterBoundary const&
61  outer,int Niter,double speed,int mode,bool Rmin)
62 {
63  RunLoadBalance(tproc,points,outer,Niter,0,speed,mode,Rmin);
64 }
65 
66 void SetLoad(Tessellation &tproc,vector<Vector2D> &points,OuterBoundary const&
67  outer,double TargetLoad,double speed,int mode,bool Rmin)
68 {
69  RunLoadBalance(tproc,points,outer,1000,TargetLoad,speed,mode,Rmin);
70 }
71 
72 #endif
Voronoi tessellation class.
Definition: VoronoiMesh.hpp:29
Abstract class for tessellation.
void SetLoad(Tessellation &tproc, vector< Vector2D > &points, OuterBoundary const &outer, int Niter=100, double speed=0.04, int mode=2, bool Rmin=false)
Corrects the load between processors based on number of cells per processor.
Definition: SetLoad.cpp:60
virtual int GetPointNo(void) const =0
Get Total number of mesh generating points.
vector< Vector2D > UpdateMPIPoints(Tessellation const &vproc, int rank, vector< Vector2D > &points, OuterBoundary const *obc, vector< size_t > &selfindex, vector< int > &sentproc, vector< vector< int > > &sentpoints)
Communicate position of mesh generating points between processes.
Function for setting the load balance based on equal cells per processor.
A load balancing scheme aiming for the same number of points in each process.
void SetPointNo(int N)
Set the number of points.
vector< Vector2D > & GetMeshPoints(void)
Returns a reference to the point vector.
Abstract class for geometric boundary conditions for the tessellation.