SpinParser  1.0
LoadManager.hpp
Go to the documentation of this file.
1 
9 #pragma once
10 #include <vector>
11 #include <functional>
12 #include <thread>
13 #include <mutex>
14 #include <boost/date_time.hpp>
15 #include "lib/Log.hpp"
16 #include "lib/Exception.hpp"
17 
18 #ifndef DISABLE_MPI
19 #include "mpi.h"
20 #endif
21 #ifndef DISABLE_OMP
22 #include "omp.h"
23 #endif
24 
25 #define HMP_CHUNK_PROPERTY_STACK 0
26 #define HMP_CHUNK_PROPERTY_BEGIN 1
27 #define HMP_CHUNK_PROPERTY_END 2
28 
29 #ifndef DISABLE_MPI
30 #define HMP_MPI_ENABLED
31 #define HMP_ENABLE_IF_MPI(X) X
32 #define HMP_DISABLE_IF_MPI(X)
33 #define HMP_APPEND_IF_MPI(X) ,X
34 #else
35 #define HMP_ENABLE_IF_MPI(X)
36 #define HMP_DISABLE_IF_MPI(X) X
37 #define HMP_APPEND_IF_MPI(X)
38 #endif
39 
40 namespace HMP
41 {
42  class LoadManager;
43  typedef int StackIdentifier;
44 
67  {
68  protected:
77  {
81  enum struct StackType
82  {
89  Explicit,
96  Implicit,
104  Slave,
110  Passive
111  };
112 
116  enum struct MessageTag
117  {
118  Chunk,
119  ChunkReturn,
120  };
121 
125  virtual ~DataStackBase() {};
126 
132  virtual void applyCalculator(const int n) {};
133 
134  #ifdef HMP_MPI_ENABLED
135 
143  virtual void send(const int offset, const int count, const int serverRank, const MPI_Comm communicator) const {};
144 
154  virtual void receive(const int offset, const int count, const int rank, const MPI_Comm communicator, MPI_Request &request) {};
155 
162  virtual void broadcast(const int serverRank, const MPI_Comm communicator) {};
163  #endif
164 
165  StackType type;
167  int size;
172  };
173 
179  template <class StackT> struct DataStack : public DataStackBase
180  {
181  friend class LoadManager;
182  protected:
192 
198  void applyCalculator(const int i) override
199  {
201  else if (type == StackType::Explicit) data[i] = explicitCalculator(i);
202  }
203 
204  #ifdef HMP_MPI_ENABLED
205 
213  void send(const int offset, const int count, const int serverRank, const MPI_Comm communicator) const override
214  {
215  MPI_Send(static_cast<void *>(data + typeMultiplicity * offset), typeMultiplicity *count * sizeof(StackT), MPI_BYTE, serverRank, static_cast<int>(MessageTag::ChunkReturn), communicator);
216  }
217 
227  void receive(const int offset, const int count, const int rank, const MPI_Comm communicator, MPI_Request &request) override
228  {
229  MPI_Irecv(static_cast<void *>(data + typeMultiplicity * offset), typeMultiplicity *count * sizeof(StackT), MPI_BYTE, rank, static_cast<int>(MessageTag::ChunkReturn), communicator, &request);
230  }
231 
238  void broadcast(const int serverRank, const MPI_Comm communicator) override
239  {
240  MPI_Bcast(static_cast<void *>(data), typeMultiplicity * size * sizeof(StackT), MPI_BYTE, serverRank, communicator);
241  }
242  #endif
243 
244  std::function<StackT(int)> explicitCalculator;
245  std::function<void(int)> implicitCalculator;
246  StackT *data;
247  };
248 
252  struct Chunk
253  {
254  public:
257  {
261  };
262 
270  Chunk(const int stackId, const int begin, const int end)
271  {
275  }
276 
278  bool isVoid() const
279  {
280  return (properties[HMP_CHUNK_PROPERTY_STACK] == -1);
281  }
282 
284  bool operator==(const Chunk &rhs) const
285  {
287  }
288 
290  bool operator!=(const Chunk &rhs) const
291  {
292  return !this->operator==(rhs);
293  }
294 
295  int properties[3];
296  };
297 
298  public:
300  virtual ~LoadManager()
301  {
302  HMP_ENABLE_IF_MPI(MPI_Comm_free(&_communicator));
303 
304  while (_stacks.size() > 0)
305  {
306  delete _stacks.back();
307  _stacks.pop_back();
308  }
309  }
310 
326  template <class StackT> StackIdentifier addMasterStackExplicit(StackT *const data, const int size, const std::function<StackT(int)> &calculator, const int recommendedChunkSizeMultiple = 1, const int recommendedChunksPerRank = 10, const bool autoBroadcast = false)
327  {
330  ds->master = -1;
331  ds->size = size;
332  ds->typeMultiplicity = 1;
333  ds->recommendedChunkSizeMultiple = recommendedChunkSizeMultiple;
334  ds->recommendedChunksPerRank = recommendedChunksPerRank;
335  ds->autoBroadcast = autoBroadcast;
336  ds->explicitCalculator = calculator;
337  ds->data = data;
338  return _registerStack(ds);
339  }
340 
357  template <class StackT> StackIdentifier addMasterStackImplicit(StackT *const data, const int size, const std::function<void(int)> &calculator, const int typeMultiplicity = 1, const int recommendedChunkSizeMultiple = 1, const int recommendedChunksPerRank = 10, const bool autoBroadcast = false)
358  {
361  ds->master = -1;
362  ds->size = size;
363  ds->typeMultiplicity = typeMultiplicity;
364  ds->recommendedChunkSizeMultiple = recommendedChunkSizeMultiple;
365  ds->recommendedChunksPerRank = recommendedChunksPerRank;
366  ds->autoBroadcast = autoBroadcast;
367  ds->implicitCalculator = calculator;
368  ds->data = data;
369  return _registerStack(ds);
370  }
371 
385  template <class StackT> StackIdentifier addSlaveStack(StackT *const data, const int size, const StackIdentifier master, const int typeMultiplicity = 1)
386  {
389  ds->master = master;
390  ds->size = size;
391  ds->typeMultiplicity = typeMultiplicity;
392  ds->autoBroadcast = false;
393  ds->data = data;
394  return _registerStack(ds);
395  }
396 
408  template <class StackT> StackIdentifier addPassiveStack(StackT *const data, const int size)
409  {
412  ds->master = -1;
413  ds->size = size;
414  ds->typeMultiplicity = 1;
415  ds->autoBroadcast = false;
416  ds->data = data;
417  return _registerStack(ds);
418  }
419 
426  virtual void calculate(const StackIdentifier *stackIds, const int size) = 0;
427 
433  void calculate(const std::initializer_list<StackIdentifier> &stackIds)
434  {
435  calculate(stackIds.begin(), int(stackIds.size()));
436  }
437 
443  void calculate(const StackIdentifier stackId)
444  {
445  calculate({ stackId });
446  }
447 
452  {
453  std::vector<StackIdentifier> all(_stacks.size());
454  for (StackIdentifier i = 0; i < StackIdentifier(_stacks.size()); ++i) all[i] = i;
455  calculate(all.data(), int(all.size()));
456  }
457 
464  void broadcast(const StackIdentifier *stackIds, const int size)
465  {
466  #ifdef HMP_MPI_ENABLED
467  for (int i = 0; i < size; ++i)
468  {
469  for (StackIdentifier s = 0; s < StackIdentifier(_stacks.size()); ++s)
470  {
471  if (s == stackIds[i] || _stacks[s]->master == stackIds[i]) _stacks[s]->broadcast(_serverRank, _communicator);
472  }
473  }
474  #endif
475  }
476 
482  void broadcast(const std::initializer_list<StackIdentifier> &stackIds)
483  {
484  broadcast(stackIds.begin(), int(stackIds.size()));
485  }
486 
492  void broadcast(const StackIdentifier stackId)
493  {
494  broadcast({ stackId });
495  }
496 
501  {
502  std::vector<StackIdentifier> all(_stacks.size());
503  for (StackIdentifier i = 0; i < StackIdentifier(_stacks.size()); ++i) all[i] = i;
504  broadcast(all.data(), int(all.size()));
505  }
506 
510  virtual void printRuntimeStatistics() const {}
511 
512  protected:
519  LoadManager(const int serverRank HMP_APPEND_IF_MPI(const MPI_Comm communicator))
520  {
521  HMP_ENABLE_IF_MPI(MPI_Comm_dup(communicator, &_communicator));
522  HMP_ENABLE_IF_MPI(MPI_Comm_size(_communicator, &_commSize));
524  HMP_ENABLE_IF_MPI(MPI_Comm_rank(_communicator, &_rank));
526 
527  if (serverRank >= _commSize) throw Exception(Exception::Type::MpiError, "Server rank must not exceed communicator size.");
528  else _serverRank = serverRank;
529  }
530 
538  {
539  //If stack is a slave stack, check that the type multiplicity matches the type multiplicity of the master
540  if (stack->type == DataStackBase::StackType::Slave && stack->typeMultiplicity != _stacks[stack->master]->typeMultiplicity) throw Exception(Exception::Type::ArgumentError, "Stack type multiplicity of slave stack does not match type multiplicity of associated master stack.");
541 
542  _stacks.push_back(stack);
543  return StackIdentifier(_stacks.size() - 1);
544  }
545 
551  void _calculateChunk(const Chunk &chunk)
552  {
553  #ifndef DISABLE_OMP
554  #pragma omp parallel for schedule(guided)
555  #endif
556  for (int i = chunk.properties[HMP_CHUNK_PROPERTY_BEGIN]; i < chunk.properties[HMP_CHUNK_PROPERTY_END]; ++i) _stacks[chunk.properties[HMP_CHUNK_PROPERTY_STACK]]->applyCalculator(i);
557  }
558 
559  std::vector<DataStackBase *> _stacks;
561  int _rank;
562  int _commSize;
564  };
565 
572  {
573  friend LoadManager *newLoadManager(const int serverRank HMP_APPEND_IF_MPI(const MPI_Comm communicator));
574  public:
581  virtual void calculate(const StackIdentifier *stackIds, const int size) override
582  {
583  //reset calculation runtime statistics
584  for (int i = 0; i < _commSize; ++i)
585  {
586  for (StackIdentifier s = 0; s < StackIdentifier(_stacks.size()); ++s)
587  {
588  _currentCalculationComputeTimeBuffer[i * _stacks.size() + s] = 0.0f;
589  _currentCalculationWorkDone[i][s] = 0;
590  _currentCalculationTime[i][s] = 0.0f;
591 
592  }
593  }
594  boost::posix_time::ptime tic = boost::posix_time::microsec_clock::local_time();
595 
596  //init chunk spawner
597  _initChunkSpawner(stackIds, size);
598 
599  //spawn local worker
600  std::thread *t = new std::thread([&]() { this->_runLocalClient(); });
601 
602  #ifdef HMP_MPI_ENABLED
603  //issue initial chunks
604  for (int i = 0; i < _commSize; ++i) if (i != _serverRank) _issueChunk(i);
605 
606  //collect results and issue consecutive chunks
607  int receiveFlag;
608  MPI_Status receiveStatus;
609  for (;;)
610  {
611  begin:
612  for (int rank = 0; rank < _commSize; ++rank)
613  {
614  if (rank == _serverRank || _pendingRequests[rank] == MPI_REQUEST_NULL) continue;
615  MPI_Test(_pendingRequests + rank, &receiveFlag, &receiveStatus);
616  if (receiveFlag == 1)
617  {
618  _despawnChunk(rank);
619  _issueChunk(receiveStatus.MPI_SOURCE);
620  }
621  }
622  for (int rank = 0; rank < _commSize; ++rank) if (_pendingRequests[rank] != MPI_REQUEST_NULL) goto begin;
623  break;
624  }
625  #endif
626 
627  //join local worker
628  t->join();
629 
630  //broadcast result
631  for (int s = 0; s < size; ++s)
632  {
633  if (_stacks[stackIds[s]]->autoBroadcast) broadcast(stackIds[s]);
634  }
635 
636  //gather runtime statistics
637  HMP_ENABLE_IF_MPI(MPI_Gather(MPI_IN_PLACE, int(_stacks.size()), MPI_FLOAT, _currentCalculationComputeTimeBuffer.data(), int(_stacks.size()), MPI_FLOAT, _serverRank, _communicator));
638 
639  for (int i = 0; i < _commSize; ++i)
640  {
641  for (StackIdentifier s = 0; s < StackIdentifier(_stacks.size()); ++s)
642  {
644  }
645  }
646 
647  boost::posix_time::ptime toc = boost::posix_time::microsec_clock::local_time();
648  _totalCalculationTime += float((toc - tic).total_milliseconds());
649  }
650 
654  void printRuntimeStatistics() const override
655  {
656  Log::log << Log::LogLevel::Debug << "LoadManager total time spent in calculate() calls is " << _totalCalculationTime << "ms" << Log::endl;
657  for (int i = 0; i < _commSize; ++i)
658  {
659  float timePerRank = 0.0;
660  for (StackIdentifier j = 0; j < StackIdentifier(_stacks.size()); ++j) timePerRank += _totalComputeTime[i][j];
661  Log::log << Log::LogLevel::Debug << "LoadManager (rank " << i << ") active computing time was " << std::setiosflags(std::ios::fixed) << std::setprecision(0) << timePerRank << "ms (Efficiency: " << std::setprecision(1) << 100.0f * timePerRank / _totalCalculationTime << "%)" << Log::endl;
662 
663  for (StackIdentifier j = 0; j < StackIdentifier(_stacks.size()); ++j) Log::log << Log::LogLevel::Debug << "\t active computing time on stack " << j << " was " << _totalComputeTime[i][j] << "ms" << Log::endl;
664  }
665  }
666 
667  protected:
674  LoadManagerMaster(const int serverRank HMP_APPEND_IF_MPI(const MPI_Comm communicator)) : LoadManager(serverRank HMP_APPEND_IF_MPI(communicator))
675  {
676  _totalCalculationTime = 0.0f;
677  _totalComputeTime = new std::vector<float>[_commSize];
678  _currentCalculationWorkDone = new std::vector<int>[_commSize];
679  _currentCalculationTime = new std::vector<float>[_commSize];
680  _currentCalculationChunkSpawntime = new boost::posix_time::ptime[_commSize];
682 
683  #ifdef HMP_MPI_ENABLED
684  _pendingRequests = new MPI_Request[_commSize];
685  for (int i = 0; i < _commSize; ++i) _pendingRequests[i] = MPI_REQUEST_NULL;
686  #endif
687  }
688 
693  {
694  delete[] _totalComputeTime;
696  delete[] _currentCalculationTime;
699 
701  }
702 
710  {
711  StackIdentifier identifier = LoadManager::_registerStack(stack);
712 
713  for (int i = 0; i < _commSize; ++i)
714  {
715  _totalComputeTime[i].push_back(0.0f);
716  _currentCalculationWorkDone[i].push_back(0);
717  _currentCalculationTime[i].push_back(0.0f);
718  }
719 
721  _currentCalculationStackMask.resize(_stacks.size());
723 
724  return identifier;
725  }
726 
731  {
732  for (;;)
733  {
734  //get chunk
736  if (c.isVoid()) break;
737 
738  boost::posix_time::ptime tic = boost::posix_time::microsec_clock::local_time();
739  _calculateChunk(c);
740  boost::posix_time::ptime toc = boost::posix_time::microsec_clock::local_time();
741  _currentCalculationComputeTimeBuffer[c.properties[HMP_CHUNK_PROPERTY_STACK]] += (toc - tic).total_milliseconds();
742 
744  }
745  }
746 
753  void _initChunkSpawner(const StackIdentifier *stackIds, const int size)
754  {
755  //reset calculation statistics
756  for (int i = 0; i < _commSize; ++i)
757  {
758  for (StackIdentifier s = 0; s < StackIdentifier(_stacks.size()); ++s)
759  {
760  _currentCalculationWorkDone[i][s] = 0;
761  _currentCalculationTime[i][s] = 0.0f;
762  }
763  }
764  for (StackIdentifier s = 0; s < StackIdentifier(_stacks.size()); ++s)
765  {
766  _currentCalculationStackMask[s] = false;
768  }
769 
770  //enable selected stacks
771  for (int i = 0; i < size; ++i) _currentCalculationStackMask[stackIds[i]] = true;
772  }
773 
779  void _initChunkSpawner(const std::initializer_list<StackIdentifier> &stackIds)
780  {
781  _initChunkSpawner(stackIds.begin(), int(stackIds.size()));
782  }
783 
796  Chunk _spawnChunk(const int rank)
797  {
798  std::lock_guard<std::mutex> lock(_currentCalculationChunkSpawnerLock);
799 
800  Chunk c;
801  for (StackIdentifier s = 0; s < StackIdentifier(_stacks.size()); ++s)
802  {
803  //skip inactive and completed stacks
804  if (!_currentCalculationStackMask[s] || _currentCalculationStackProgress[s] >= _stacks[s]->size) continue;
805 
806  //determine max chunk size
807  int maximumWorkShare = int(_stacks[s]->size / (_stacks[s]->recommendedChunksPerRank * _commSize));
808  maximumWorkShare = (int(maximumWorkShare / _stacks[s]->recommendedChunkSizeMultiple) + 1) * _stacks[s]->recommendedChunkSizeMultiple;
809  if (maximumWorkShare < 1) maximumWorkShare = 1;
810 
811  //determine dynamic chunk size according to compute power
812  float myComputePower = _currentCalculationWorkDone[rank][s] / _currentCalculationTime[rank][s];
813  float totalComputePower = 0;
814  for (int i = 0; i < _commSize; ++i)
815  for (StackIdentifier j = 0; j < StackIdentifier(_stacks.size()); ++j) totalComputePower += _currentCalculationWorkDone[i][j] / _currentCalculationTime[i][j];
816 
817  int newCurrentCalculationProgress;
818  if (std::isfinite(myComputePower) && std::isfinite(totalComputePower))
819  {
820  //factor in amount of remaining work
821  int remainingWork = int(_stacks[s]->size - _currentCalculationStackProgress[s]);
822  int myWorkShare = int((myComputePower / totalComputePower) * remainingWork);
823  myWorkShare = (int(myWorkShare / _stacks[s]->recommendedChunkSizeMultiple) + 1) * _stacks[s]->recommendedChunkSizeMultiple;
824 
825  //clip to min/max chunk size
826  int minumumWorkTime = 100;
827  int minimumWorkShare = int(myComputePower * minumumWorkTime);
828  if (minimumWorkShare < 1) minimumWorkShare = 1;
829  if (myWorkShare < minimumWorkShare) myWorkShare = minimumWorkShare;
830  if (myWorkShare > maximumWorkShare) myWorkShare = maximumWorkShare;
831  newCurrentCalculationProgress = _currentCalculationStackProgress[s] + myWorkShare;
832  }
833  else newCurrentCalculationProgress = _currentCalculationStackProgress[s] + maximumWorkShare;
834 
835  //clip chunk to stack size
836  if (newCurrentCalculationProgress >= _stacks[s]->size) newCurrentCalculationProgress = _stacks[s]->size;
837 
838  //write chunk parameters
839  _currentCalculationWorkDone[rank][s] += newCurrentCalculationProgress - _currentCalculationStackProgress[s];
842  c.properties[HMP_CHUNK_PROPERTY_END] = decltype(c.properties[HMP_CHUNK_PROPERTY_END])(newCurrentCalculationProgress);
843  Log::log << Log::LogLevel::Debug << "LoadManager spawned chunk (stack " << s << ", from " << _currentCalculationStackProgress[s] << ", to " << newCurrentCalculationProgress << ", rank " << rank << ")" << Log::endl;
844 
845  //return chunk
846  _currentCalculationStackProgress[s] = newCurrentCalculationProgress;
847  break;
848  }
849  _currentCalculationChunkSpawntime[rank] = boost::posix_time::microsec_clock::local_time();
851  return c;
852  }
853 
859  void _issueChunk(const int rank)
860  {
861  #ifdef HMP_MPI_ENABLED
862  Chunk c = _spawnChunk(rank);
863  MPI_Send(&c.properties, 3, MPI_INT, rank, static_cast<int>(DataStackBase::MessageTag::Chunk), _communicator);
864 
865  if (!c.isVoid())
866  {
867  for (StackIdentifier i = 0; i < StackIdentifier(_stacks.size()); ++i)
868  {
869  //we may expect to receive data from additional slave stacks; due to MPI's non-overtaking policy, it is sufficient to only store the most recent request per rank
871  }
872  }
873  else _pendingRequests[rank] = MPI_REQUEST_NULL;
874  #endif
875  }
876 
882  void _despawnChunk(const int rank)
883  {
884  float chunktime = float((boost::posix_time::microsec_clock::local_time() - _currentCalculationChunkSpawntime[rank]).total_milliseconds());
885  std::lock_guard<std::mutex> lock(_currentCalculationChunkSpawnerLock);
887  }
888 
890  std::vector<float> *_totalComputeTime;
891  std::vector<int> *_currentCalculationWorkDone;
892  std::vector<float> *_currentCalculationTime;
893  boost::posix_time::ptime *_currentCalculationChunkSpawntime;
896  std::vector<bool> _currentCalculationStackMask;
899 
901  };
902 
909  {
910  friend LoadManager *newLoadManager(const int serverRank HMP_APPEND_IF_MPI(const MPI_Comm communicator));
911  public:
918  virtual void calculate(const StackIdentifier *stackIds, const int size) override
919  {
920  #ifdef HMP_MPI_ENABLED
921  memset(_currentCalculationComputeTimeBuffer.data(), 0, _currentCalculationComputeTimeBuffer.size() * sizeof(float));
922 
923  Chunk c;
924  for (;;)
925  {
926  //receive chunk
927  _waitChunk(c);
928  if (c.isVoid()) break;
929 
930  //compute chunk
931  boost::posix_time::ptime tic = boost::posix_time::microsec_clock::local_time();
932  _calculateChunk(c);
933  boost::posix_time::ptime toc = boost::posix_time::microsec_clock::local_time();
934  _currentCalculationComputeTimeBuffer[c.properties[HMP_CHUNK_PROPERTY_STACK]] += (toc - tic).total_milliseconds();
935 
936  //return chunk
937  _returnChunk(c);
938  }
939 
940  //broadcast result
941  for (int s = 0; s < size; ++s)
942  {
943  if (_stacks[stackIds[s]]->autoBroadcast) broadcast(stackIds[s]);
944  }
945 
946  //gather compute time statistics
947  MPI_Gather(_currentCalculationComputeTimeBuffer.data(), int(_currentCalculationComputeTimeBuffer.size()), MPI_FLOAT, nullptr, int(_currentCalculationComputeTimeBuffer.size()), MPI_FLOAT, _serverRank, _communicator);
948  #endif
949  }
950 
951  protected:
958  LoadManagerSlave(const int serverRank HMP_APPEND_IF_MPI(const MPI_Comm communicator)) : LoadManager(serverRank HMP_APPEND_IF_MPI(communicator)) {}
959 
967  {
968  StackIdentifier identifier = LoadManager::_registerStack(stack);
969 
971 
972  return identifier;
973  }
974 
980  void _waitChunk(Chunk &chunk) const
981  {
982  #ifdef HMP_MPI_ENABLED
983  MPI_Status status;
984  MPI_Recv(&chunk.properties, 3, MPI_INT, _serverRank, static_cast<int>(DataStackBase::MessageTag::Chunk), _communicator, &status);
985  #endif
986  }
987 
993  void _returnChunk(const Chunk &chunk) const
994  {
995  #ifdef HMP_MPI_ENABLED
996  for (int i = 0; i < int(_stacks.size()); ++i)
997  {
999  }
1000  #endif
1001  }
1002 
1004  };
1005 
1014  inline LoadManager *newLoadManager(const int serverRank = 0 HMP_APPEND_IF_MPI(const MPI_Comm communicator = MPI_COMM_WORLD))
1015  {
1016  int myRank;
1017  HMP_ENABLE_IF_MPI(MPI_Comm_rank(communicator, &myRank));
1018  HMP_DISABLE_IF_MPI(myRank = 0);
1019 
1020  if (myRank == serverRank) return new LoadManagerMaster(serverRank HMP_APPEND_IF_MPI(communicator));
1021  else return new LoadManagerSlave(serverRank HMP_APPEND_IF_MPI(communicator));
1022  };
1023 
1024 } //namespace HMP
1025 
1026 #undef HMP_CHUNK_PROPERTY_STACK
1027 #undef HMP_CHUNK_PROPERTY_BEGIN
1028 #undef HMP_CHUNK_PROPERTY_END
1029 
1030 #undef HMP_MPI_ENABLED
1031 #undef HMP_ENABLE_IF_MPI
1032 #undef HMP_DISABLE_IF_MPI
1033 #undef HMP_APPEND_IF_MPI
HMP::LoadManager::DataStackBase::StackType::Passive
@ Passive
Passive DataStack.
HMP::LoadManager::DataStackBase::recommendedChunkSizeMultiple
int recommendedChunkSizeMultiple
When breaking the data stack down into smaller work chunks, attempt to form chunks whose size is a mu...
Definition: LoadManager.hpp:169
HMP::LoadManager::Chunk
Definition of a small chunk of work. Used to communicate workload between different LoadManager insta...
Definition: LoadManager.hpp:252
HMP::LoadManagerMaster
LoadManager master implementation, which is responsible for distributing work and synchronizing data.
Definition: LoadManager.hpp:571
HMP::LoadManager::calculateAll
void calculateAll()
Calculate all stacks.
Definition: LoadManager.hpp:451
HMP::LoadManagerMaster::_totalComputeTime
std::vector< float > * _totalComputeTime
_totalComputeTime[rank][stack] is the accumulated time in milliseconds which MPI rank rank spent comp...
Definition: LoadManager.hpp:890
HMP::LoadManager::printRuntimeStatistics
virtual void printRuntimeStatistics() const
Virtual implementation to print runtime statistics, including information on the efficiency of LoadMa...
Definition: LoadManager.hpp:510
HMP::LoadManager::_communicator
MPI_Comm _communicator
MPI communicator to operate on.
Definition: LoadManager.hpp:563
HMP::LoadManager::DataStack::data
StackT * data
Internal data array.
Definition: LoadManager.hpp:246
HMP::LoadManager::DataStackBase::~DataStackBase
virtual ~DataStackBase()
Virtual destructor.
Definition: LoadManager.hpp:125
HMP::LoadManager::addSlaveStack
StackIdentifier addSlaveStack(StackT *const data, const int size, const StackIdentifier master, const int typeMultiplicity=1)
Create a slave data stack and attach it to the LoadManager.
Definition: LoadManager.hpp:385
HMP::LoadManager::DataStack::receive
void receive(const int offset, const int count, const int rank, const MPI_Comm communicator, MPI_Request &request) override
Asynchronously receive a data block from a different MPI rank.
Definition: LoadManager.hpp:227
HMP::LoadManager::DataStackBase::recommendedChunksPerRank
int recommendedChunksPerRank
When breaking the data stack down into smaller work chunks, attempt to form approximately the specifi...
Definition: LoadManager.hpp:170
HMP::LoadManager::DataStack::explicitCalculator
std::function< StackT(int)> explicitCalculator
Explicit calculator. Only relevant if DataStackBase::type is set to StackType::Explicit.
Definition: LoadManager.hpp:244
HMP_ENABLE_IF_MPI
#define HMP_ENABLE_IF_MPI(X)
Print argument if MPI parallelization is enabled.
Definition: LoadManager.hpp:31
HMP::LoadManager::_registerStack
virtual StackIdentifier _registerStack(DataStackBase *stack)
Register a DataStackBase with the LoadManager. By registering the stack, the LoadManager assumes resp...
Definition: LoadManager.hpp:537
HMP_CHUNK_PROPERTY_STACK
#define HMP_CHUNK_PROPERTY_STACK
Memory offset of the stack id in the chunk properties.
Definition: LoadManager.hpp:25
HMP::StackIdentifier
int StackIdentifier
DataStack identifier.
Definition: LoadManager.hpp:42
HMP::LoadManagerMaster::LoadManagerMaster
LoadManagerMaster(const int serverRank, const MPI_Comm communicator)
Construct a new LoadManagerMaster object.
Definition: LoadManager.hpp:674
HMP::LoadManagerMaster::_spawnChunk
Chunk _spawnChunk(const int rank)
The chunk spawner routine, which is responsible for breaking down the overall workload of calculator ...
Definition: LoadManager.hpp:796
HMP::LoadManager::DataStack::DataStack
DataStack()
Construct a new Data Stack object.
Definition: LoadManager.hpp:191
HMP::LoadManager::DataStackBase::size
int size
Number of elements (or tuples of elements if DataStackBase::typeMultiplicity is greater than one) sto...
Definition: LoadManager.hpp:167
HMP::LoadManagerMaster::printRuntimeStatistics
void printRuntimeStatistics() const override
Print runtime statistics, including information on the efficiency of LoadManager instances running on...
Definition: LoadManager.hpp:654
HMP::LoadManager::DataStackBase::StackType
StackType
List to identify the different traits of DataStack.
Definition: LoadManager.hpp:81
HMP::LoadManager::broadcast
void broadcast(const StackIdentifier stackId)
Broadcast a single stack.
Definition: LoadManager.hpp:492
HMP::LoadManager::_stacks
std::vector< DataStackBase * > _stacks
List of all registered stacks.
Definition: LoadManager.hpp:559
HMP::LoadManagerMaster::~LoadManagerMaster
~LoadManagerMaster()
Destroy the LoadManagerMaster object.
Definition: LoadManager.hpp:692
HMP::LoadManager::DataStackBase::StackType::Slave
@ Slave
Slave DataStack.
HMP::newLoadManager
LoadManager * newLoadManager(const int serverRank=0, const MPI_Comm communicator=MPI_COMM_WORLD)
Create a new LoadManager instance. If the current MPI rank is the designated master rank,...
Definition: LoadManager.hpp:1014
HMP::LoadManagerMaster::_currentCalculationChunkSpawntime
boost::posix_time::ptime * _currentCalculationChunkSpawntime
_currentCalculationChunkSpawntime[rank] specifies time at which the most recent chunk has been issued...
Definition: LoadManager.hpp:893
HMP::LoadManager::DataStackBase::MessageTag::ChunkReturn
@ ChunkReturn
MPI message contains the result of a chunk computation.
HMP::LoadManager::_rank
int _rank
MPI rank of the current LoadManager instance.
Definition: LoadManager.hpp:561
HMP::LoadManager::DataStack::broadcast
void broadcast(const int serverRank, const MPI_Comm communicator) override
Broadcast a data block from the server rank to all other MPI ranks.
Definition: LoadManager.hpp:238
HMP::LoadManagerSlave::_returnChunk
void _returnChunk(const Chunk &chunk) const
Return the result of the workload defined by a specific chunk.
Definition: LoadManager.hpp:993
HMP::LoadManagerMaster::_pendingRequests
MPI_Request * _pendingRequests
MPI request objects associated with the return values for workload chunks that have been issued.
Definition: LoadManager.hpp:900
HMP::LoadManager::calculate
virtual void calculate(const StackIdentifier *stackIds, const int size)=0
Calculate a list of stacks, where the stack identifiers are provided in list form.
HMP::LoadManager::calculate
void calculate(const std::initializer_list< StackIdentifier > &stackIds)
Calculate a list of stacks, where the stack identifiers are provided in initializer list form.
Definition: LoadManager.hpp:433
HMP::LoadManagerMaster::_currentCalculationStackProgress
std::vector< int > _currentCalculationStackProgress
_currentCalculationStackProgress[stack] specifies the current progress (pointer to the next unissued ...
Definition: LoadManager.hpp:897
HMP::LoadManager::DataStackBase::autoBroadcast
bool autoBroadcast
If set to true, modifications to the stack's data that are a consequence of the onvication of calcula...
Definition: LoadManager.hpp:171
Exception::Type::MpiError
@ MpiError
MPI error, raised when MPI communication fails.
HMP::LoadManagerMaster::_currentCalculationWorkDone
std::vector< int > * _currentCalculationWorkDone
_currentCalculationWorkDone[rank][stack] is the number of calculations which have been performed by M...
Definition: LoadManager.hpp:891
HMP::LoadManagerSlave::calculate
virtual void calculate(const StackIdentifier *stackIds, const int size) override
Calculate a list of stacks, where the stack identifiers are provided in list form.
Definition: LoadManager.hpp:918
HMP::LoadManagerMaster::_currentCalculationStackMask
std::vector< bool > _currentCalculationStackMask
_currentCalculationStackMask[stack] specifies whether stack should be computed in the current calcula...
Definition: LoadManager.hpp:896
HMP::LoadManager::calculate
void calculate(const StackIdentifier stackId)
Calculate a single stack.
Definition: LoadManager.hpp:443
HMP::LoadManager::_calculateChunk
void _calculateChunk(const Chunk &chunk)
Calculate the workload defined by a specific chunk.
Definition: LoadManager.hpp:551
Exception
Descriptor object for exceptions.
Definition: Exception.hpp:17
HMP::LoadManager::Chunk::properties
int properties[3]
Workload specification. First value describes the stack id, second value describes the first element ...
Definition: LoadManager.hpp:295
HMP_CHUNK_PROPERTY_BEGIN
#define HMP_CHUNK_PROPERTY_BEGIN
Memory offset of the workload begin in the chunk properties.
Definition: LoadManager.hpp:26
HMP::LoadManagerMaster::_currentCalculationChunkSpawned
Chunk * _currentCalculationChunkSpawned
_currentCalculationChunkSpawned[rank] stores the most recent chunk generated for MPI rank rank in the...
Definition: LoadManager.hpp:894
HMP::LoadManager::Chunk::operator!=
bool operator!=(const Chunk &rhs) const
Negative comparison operator.
Definition: LoadManager.hpp:290
HMP::LoadManager::addMasterStackExplicit
StackIdentifier addMasterStackExplicit(StackT *const data, const int size, const std::function< StackT(int)> &calculator, const int recommendedChunkSizeMultiple=1, const int recommendedChunksPerRank=10, const bool autoBroadcast=false)
Create an explicit data stack and attach it to the LoadManager.
Definition: LoadManager.hpp:326
HMP::LoadManagerMaster::_issueChunk
void _issueChunk(const int rank)
Send a workload chunk to a specified MPI rank.
Definition: LoadManager.hpp:859
HMP::LoadManager::DataStackBase::StackType::Implicit
@ Implicit
Implicit DataStack.
HMP::LoadManagerMaster::calculate
virtual void calculate(const StackIdentifier *stackIds, const int size) override
Calculate a list of stacks, where the stack identifiers are provided in list form.
Definition: LoadManager.hpp:581
HMP::LoadManager::DataStackBase::broadcast
virtual void broadcast(const int serverRank, const MPI_Comm communicator)
Virtual function to broadcast a data block from the server rank to all other MPI ranks.
Definition: LoadManager.hpp:162
HMP::LoadManager::broadcast
void broadcast(const StackIdentifier *stackIds, const int size)
Broadcast a list of stacks, where the stack identifiers are provided in list form.
Definition: LoadManager.hpp:464
HMP::LoadManager::_commSize
int _commSize
MPI communicator size.
Definition: LoadManager.hpp:562
HMP::LoadManagerSlave::LoadManagerSlave
LoadManagerSlave(const int serverRank, const MPI_Comm communicator)
Construct a new LoadManagerSlave object.
Definition: LoadManager.hpp:958
HMP::LoadManager::LoadManager
LoadManager(const int serverRank, const MPI_Comm communicator)
Construct a new Load Manager object.
Definition: LoadManager.hpp:519
HMP::LoadManagerSlave
LoadManager slave implementation, responsible for receiving and executing workload chunks form a Load...
Definition: LoadManager.hpp:908
HMP::LoadManagerSlave::newLoadManager
friend LoadManager * newLoadManager(const int serverRank, const MPI_Comm communicator)
Create a new LoadManager instance. If the current MPI rank is the designated master rank,...
Definition: LoadManager.hpp:1014
HMP::LoadManagerMaster::_initChunkSpawner
void _initChunkSpawner(const StackIdentifier *stackIds, const int size)
Prepare the chunk spawner, i.e. the routine responsible for breaking down workload into smaller chuns...
Definition: LoadManager.hpp:753
HMP::LoadManager::~LoadManager
virtual ~LoadManager()
Destroy the LoadManager object.
Definition: LoadManager.hpp:300
HMP::LoadManagerSlave::_currentCalculationComputeTimeBuffer
std::vector< float > _currentCalculationComputeTimeBuffer
_currentCalculationComputeTimeBuffer[stack] is a buffer for the time in milliseconds spent on computi...
Definition: LoadManager.hpp:1003
HMP::LoadManager::Chunk::isVoid
bool isVoid() const
Check of the workload is empty.
Definition: LoadManager.hpp:278
HMP::LoadManager::Chunk::operator==
bool operator==(const Chunk &rhs) const
Comparison operator.
Definition: LoadManager.hpp:284
HMP::LoadManager::DataStackBase::receive
virtual void receive(const int offset, const int count, const int rank, const MPI_Comm communicator, MPI_Request &request)
Virtual function to asynchronously receive a data block from a different MPI rank.
Definition: LoadManager.hpp:154
HMP::LoadManager::DataStackBase
Abstract base class for DataStack types.
Definition: LoadManager.hpp:76
HMP_CHUNK_PROPERTY_END
#define HMP_CHUNK_PROPERTY_END
Memory offset of the workload end in the chunk properties.
Definition: LoadManager.hpp:27
HMP::LoadManagerMaster::_registerStack
virtual StackIdentifier _registerStack(DataStackBase *stack) override
Register a DataStackBase with the LoadManager. By registering the stack, the LoadManager assumes resp...
Definition: LoadManager.hpp:709
HMP::LoadManager::DataStackBase::type
StackType type
Specifies the trait of the stack.
Definition: LoadManager.hpp:162
Exception::Type::ArgumentError
@ ArgumentError
Argument error, raised when a function is invoked with an invalid argument.
HMP::LoadManagerMaster::_totalCalculationTime
float _totalCalculationTime
Accumulated time in milliseconds which has been spent on calculate() calls over the lifetime of the L...
Definition: LoadManager.hpp:889
HMP::LoadManagerMaster::_initChunkSpawner
void _initChunkSpawner(const std::initializer_list< StackIdentifier > &stackIds)
Prepare the chunk spawner, i.e. the routine responsible for breaking down workload into smaller chuns...
Definition: LoadManager.hpp:779
HMP::LoadManagerSlave::_registerStack
virtual StackIdentifier _registerStack(DataStackBase *stack) override
Register a DataStackBase with the LoadManager. By registering the stack, the LoadManager assumes resp...
Definition: LoadManager.hpp:966
HMP::LoadManagerMaster::_runLocalClient
void _runLocalClient()
Worker loop to run calculators locally.
Definition: LoadManager.hpp:730
HMP::LoadManager::addMasterStackImplicit
StackIdentifier addMasterStackImplicit(StackT *const data, const int size, const std::function< void(int)> &calculator, const int typeMultiplicity=1, const int recommendedChunkSizeMultiple=1, const int recommendedChunksPerRank=10, const bool autoBroadcast=false)
Create an implicit data stack and attach it to the LoadManager.
Definition: LoadManager.hpp:357
HMP::LoadManager::DataStackBase::applyCalculator
virtual void applyCalculator(const int n)
Virtual function to apply an internal calculator to the n-th value stored in the stack.
Definition: LoadManager.hpp:132
HMP::LoadManager::DataStackBase::MessageTag::Chunk
@ Chunk
MPI message contains Chunk information.
Exception.hpp
Descriptor object for exceptions.
HMP::LoadManagerMaster::_currentCalculationChunkSpawnerLock
std::mutex _currentCalculationChunkSpawnerLock
Lock to synchronize chunk spawning for remote calculations and for local worker threads.
Definition: LoadManager.hpp:898
HMP_DISABLE_IF_MPI
#define HMP_DISABLE_IF_MPI(X)
Do not print argument if MPI parallelization is enabled.
Definition: LoadManager.hpp:32
HMP::LoadManager::DataStack::send
void send(const int offset, const int count, const int serverRank, const MPI_Comm communicator) const override
Send a data block to a different MPI rank, typically the server rank.
Definition: LoadManager.hpp:213
HMP::LoadManagerMaster::_currentCalculationComputeTimeBuffer
std::vector< float > _currentCalculationComputeTimeBuffer
_currentCalculationComputeTimeBuffer[rank*_stacks.size()+stack] is a buffer for the time in milliseco...
Definition: LoadManager.hpp:895
HMP::LoadManager::DataStack::applyCalculator
void applyCalculator(const int i) override
Invoke internal calculator and write result to the data array at the specified index.
Definition: LoadManager.hpp:198
Log.hpp
Lightweight logging interface with output filtering.
HMP::LoadManager::DataStackBase::MessageTag
MessageTag
Internal message tags for MPI communication.
Definition: LoadManager.hpp:116
HMP::LoadManager::DataStack::implicitCalculator
std::function< void(int)> implicitCalculator
Implicit calculator. Only relevant if DataStackBase::type is set to StackType::Implicit.
Definition: LoadManager.hpp:245
HMP::LoadManager::DataStackBase::StackType::Explicit
@ Explicit
Explicit DataStack.
HMP::LoadManager
Common load manager interface for both, the MPI server rank and slave ranks.
Definition: LoadManager.hpp:66
HMP::LoadManagerMaster::_despawnChunk
void _despawnChunk(const int rank)
Mark a workload chunk as completed. Relevant only for runtime statistics information.
Definition: LoadManager.hpp:882
HMP::LoadManagerMaster::newLoadManager
friend LoadManager * newLoadManager(const int serverRank, const MPI_Comm communicator)
Create a new LoadManager instance. If the current MPI rank is the designated master rank,...
Definition: LoadManager.hpp:1014
HMP::LoadManager::DataStackBase::send
virtual void send(const int offset, const int count, const int serverRank, const MPI_Comm communicator) const
Virtual function to send a data block to a DataStack on a different MPI rank, typically the server ra...
Definition: LoadManager.hpp:143
HMP::LoadManager::DataStack
Concrete derivation of DataStackBase, managing an array of fundamental data types (or tuples thereof,...
Definition: LoadManager.hpp:179
HMP::LoadManager::broadcastAll
void broadcastAll()
Broadcast all stacks.
Definition: LoadManager.hpp:500
HMP::LoadManagerMaster::_currentCalculationTime
std::vector< float > * _currentCalculationTime
_currentCalculationTime[rank][stack] is the time in milliseconds spent by MPI rank rank until returni...
Definition: LoadManager.hpp:892
HMP::LoadManager::broadcast
void broadcast(const std::initializer_list< StackIdentifier > &stackIds)
Broadcast a list of stacks, where the stack identifiers are provided in initializer list form.
Definition: LoadManager.hpp:482
HMP::LoadManagerSlave::_waitChunk
void _waitChunk(Chunk &chunk) const
Wait until a workload chunk has been received from a LoadManagerMaster instance.
Definition: LoadManager.hpp:980
HMP::LoadManager::Chunk::Chunk
Chunk()
Construct an empty work chunk.
Definition: LoadManager.hpp:256
HMP::LoadManager::Chunk::Chunk
Chunk(const int stackId, const int begin, const int end)
Construct a work chunk for a specific workload.
Definition: LoadManager.hpp:270
HMP::LoadManager::_serverRank
int _serverRank
Designated MPI master rank.
Definition: LoadManager.hpp:560
HMP_APPEND_IF_MPI
#define HMP_APPEND_IF_MPI(X)
Append argument if MPI parallelizatino is enabled.
Definition: LoadManager.hpp:33
HMP::LoadManager::addPassiveStack
StackIdentifier addPassiveStack(StackT *const data, const int size)
Create a passive data stack and attach it to the LoadManager.
Definition: LoadManager.hpp:408
HMP::LoadManager::DataStackBase::typeMultiplicity
int typeMultiplicity
Multiplicity of each element. Cannot be greater than one for explicit stacks. If greater than one,...
Definition: LoadManager.hpp:168
HMP::LoadManager::DataStackBase::master
StackIdentifier master
Specifies the id of an associated stack. Only relevant for slave stacks, otherwise set to -1.
Definition: LoadManager.hpp:166