GemaCoreLib
The GeMA Core library
gmTaskManager.h
Go to the documentation of this file.
1 /************************************************************************
2 **
3 ** Copyright (C) 2014 by Carlos Augusto Teixera Mendes
4 ** All rights reserved.
5 **
6 ** This file is part of the "GeMA" software. It's use should respect
7 ** the terms in the license agreement that can be found together
8 ** with this source code.
9 ** It is provided AS IS, with NO WARRANTY OF ANY KIND,
10 ** INCLUDING THE WARRANTY OF DESIGN, MERCHANTABILITY AND FITNESS FOR
11 ** A PARTICULAR PURPOSE.
12 **
13 ************************************************************************/
14 
24 #ifndef _GEMA_TASK_MANAGER_H_
25 #define _GEMA_TASK_MANAGER_H_
26 
27 #include "gmCoreConfig.h"
28 #include "gmValueInfo.h"
29 #include "gmThreadManager.h"
30 
31 #include "gmTrace.h"
32 
33 #include <luaEnv.h>
34 
35 class GmSimulationData;
36 class GmMesh;
37 class GmCellMesh;
38 class GmCell;
39 class GmCellGroupSet;
40 
41 
44 {
45 public:
47  virtual ~GmNodeTaskDataSet() {}
48 
55  virtual int nextNode() = 0;
56 
58  virtual void reset() = 0;
59 };
60 
63 {
64 public:
66  virtual ~GmCellTaskDataSet() {}
67 
69  virtual GmCell* nextCell() = 0;
70 
72  virtual void reset() = 0;
73 };
74 
77 {
78 public:
79 
82  {
83  NPS_DEFAULT = -1,
84  NPS_SEQUENTIAL_BLOCK = 0,
85 
86  //-----------------------------
87  // No entries should be added after this line
88  // IMPORTANT: When adding entries or CHANGING enum order, please remember
89  // to also change str2NodePartitionStrategy() AND
90  // GmLuaUtils::checkAndLoadTaskOptions()
91  //-----------------------------
92  NPS_NUM_STRATEGIES
93  };
94 
97  {
98  CPS_DEFAULT = -1,
99  CPS_SEQUENTIAL_BLOCK = 0,
100 
101  //-----------------------------
102  // No entries should be added after this line
103  // IMPORTANT: When adding entries or CHANGING enum order, please remember
104  // to also change str2CellPartitionStrategy() AND
105  // GmLuaUtils::checkAndLoadTaskOptions()
106  //-----------------------------
107  CPS_NUM_STRATEGIES
108  };
109 
111 
112  static NodePartitionStrategy str2NodePartitionStrategy(QString str);
113  static CellPartitionStrategy str2CellPartitionStrategy(QString str);
114 
116  GmThreadManager* threadManager() const { return _threadManager; }
117 
119  int numTasks() const { return _threadManager->numTasks(); }
120 
121  template <class T>
122  GmThreadTaskResult runParallelNodeLoop(const GmMesh* m, GmAffectedNodes affNodes, T& taskObj,
123  int ntasks = 0, NodePartitionStrategy strategy = NPS_DEFAULT,
124  int nworkers = -1, const QVector<int>& taskAffinity = QVector<int>());
125  template <class T>
126  void execParallelNodeLoop(const GmMesh* m, GmAffectedNodes affNodes, T& taskObj,
127  int ntasks = 0, NodePartitionStrategy strategy = NPS_DEFAULT,
128  int nworkers = -1, const QVector<int>& taskAffinity = QVector<int>());
129 
130  template <class T>
131  GmThreadTaskResult runParallelCellLoop(GmCellMesh* m, bool activeOnly, T& taskObj,
132  int ntasks = 0, CellPartitionStrategy strategy = CPS_DEFAULT,
133  int nworkers = -1, const QVector<int>& taskAffinity = QVector<int>());
134 
135  template <class T>
136  GmThreadTaskResult runParallelCellLoop(GmCellGroupSet* gs, bool activeOnly, T& taskObj,
137  int ntasks = 0, CellPartitionStrategy strategy = CPS_DEFAULT,
138  int nworkers = -1, const QVector<int>& taskAffinity = QVector<int>());
139 
140  template <class T>
141  void execParallelCellLoop(GmCellMesh* m, bool activeOnly, T& taskObj,
142  int ntasks = 0, CellPartitionStrategy strategy = CPS_DEFAULT,
143  int nworkers = -1, const QVector<int>& taskAffinity = QVector<int>());
144 
145  template <class T> bool addNodeTasks(const GmMesh* m, GmAffectedNodes affNodes, T& taskObj,
146  int ntasks, NodePartitionStrategy strategy, int nworkers,
147  const QVector<int>& taskAffinity, bool* ok);
148 
149  template <class M, class T> bool addCellTasks(M* m, bool activeOnly, T& taskObj,
150  int ntasks, CellPartitionStrategy strategy, int nworkers,
151  const QVector<int>& taskAffinity, bool* ok);
152 
153 private:
154  void adjustExecParameters(int* ntasks, int* nworkers);
155 
156  QVector<GmNodeTaskDataSet*> buildNodePartition(const GmMesh* m, int n, GmAffectedNodes affNodes,
157  NodePartitionStrategy strategy) const;
158 
159  template <class M>
160  QVector<GmCellTaskDataSet*> buildCellPartition(M* m, int n, bool activeOnly,
161  CellPartitionStrategy strategy) const;
162 
164 
168 };
169 
204 template <class T>
206  int ntasks, NodePartitionStrategy strategy,
207  int nworkers, const QVector<int>& taskAffinity)
208 {
209  S_TRACE();
210 
211  // Add tasks to the manager.
212  bool ok;
213  if(!addNodeTasks<T>(m, affNodes, taskObj, ntasks, strategy, nworkers, taskAffinity, &ok))
214  return ok ? GM_THREAD_OK : GM_THREAD_ABORTED; // false + ok = no tasks
215 
216  // Run tasks!
217  return _threadManager->runTasks(nworkers);
218 }
219 
221 template <class T>
222 void GmTaskManager::execParallelNodeLoop(const GmMesh* m, GmAffectedNodes affNodes, T& taskObj,
223  int ntasks, NodePartitionStrategy strategy,
224  int nworkers, const QVector<int>& taskAffinity)
225 {
226  S_TRACE();
227 
228  // Add tasks to the manager.
229  bool ok;
230  if(!addNodeTasks<T>(m, affNodes, taskObj, ntasks, strategy, nworkers, taskAffinity, &ok))
231  {
232  if(!ok)
233  luaL_error(_threadManager->threadLuaEnv(0)->state(), "%s", luaPrintable(QObject::tr("Task execution aborted.")));
234  return;
235  }
236 
237  // Run tasks!
238  _threadManager->execTasks(nworkers);
239 }
240 
248 template <class T>
249 bool GmTaskManager::addNodeTasks(const GmMesh* m, GmAffectedNodes affNodes, T& taskObj,
250  int ntasks, NodePartitionStrategy strategy, int nworkers,
251  const QVector<int>& taskAffinity, bool* ok)
252 {
253  S_TRACE();
254 
255  *ok = true;
256 
257  // Adjust execution parameters, converting possible default options to real values
258  adjustExecParameters(&ntasks, &nworkers);
259  if(strategy == NPS_DEFAULT)
260  strategy = _defNodeStrategy;
261  assert(ntasks > 0 && nworkers >= 0 && strategy >= 0 && strategy < NPS_NUM_STRATEGIES);
262 
263  // Partition the mesh nodes
264  QVector<GmNodeTaskDataSet*> taskData = buildNodePartition(m, ntasks, affNodes, strategy);
265  if(taskData.isEmpty()) // No tasks
266  return false;
267 
268  // Check affinity vector
269  bool useAffinity = !taskAffinity.isEmpty();
270  if(useAffinity && taskData.size() != taskAffinity.size())
271  {
272  qDeleteAll(taskData);
273  gmErrorMsg(_threadManager->logger(), QObject::tr("Error creating node tasks. Affinity vector size (%1) is different from the number of tasks (%2).")
274  .arg(taskAffinity.size()).arg(taskData.size()));
275  *ok = false;
276  return false; // No way to support the affinity request
277  }
278 
279  // Add tasks to the manager. The thread manager takes ownership and destroys the tasks
280  class NodeTask: public GmThreadTask
281  {
282  public:
283  NodeTask(T& task, GmNodeTaskDataSet* dataSet) : _task(task), _dataSet(dataSet) {}
284  virtual ~NodeTask() { delete _dataSet; }
285 
286  virtual GmThreadTaskResult run(const GmThread* thread) { return _task(_dataSet, thread); }
287 
288  private:
289  T& _task;
290  GmNodeTaskDataSet* _dataSet;
291  };
292 
293  int i = 0;
294  foreach(GmNodeTaskDataSet* ds, taskData)
295  _threadManager->addTask(new NodeTask(taskObj, ds), useAffinity ? taskAffinity[i++] : 0);
296 
297  return true;
298 }
299 
334 template <class T>
336  int ntasks, CellPartitionStrategy strategy,
337  int nworkers, const QVector<int>& taskAffinity)
338 {
339  S_TRACE();
340 
341  // Add tasks to the manager.
342  bool ok;
343  if(!addCellTasks<GmCellMesh, T>(m, activeOnly, taskObj, ntasks, strategy, nworkers, taskAffinity, &ok))
344  return ok ? GM_THREAD_OK : GM_THREAD_ABORTED; // false + ok = no tasks
345 
346  // Run tasks!
347  return _threadManager->runTasks(nworkers);
348 }
349 
351 template <class T>
353  int ntasks, CellPartitionStrategy strategy,
354  int nworkers, const QVector<int>& taskAffinity)
355 {
356  S_TRACE();
357 
358  // Add tasks to the manager.
359  bool ok;
360  if(!addCellTasks<GmCellGroupSet, T>(gs, activeOnly, taskObj, ntasks, strategy, nworkers, taskAffinity, &ok))
361  return ok ? GM_THREAD_OK : GM_THREAD_ABORTED; // false + ok = no tasks
362 
363  // Run tasks!
364  return _threadManager->runTasks(nworkers);
365 }
366 
368 template <class T>
369 void GmTaskManager::execParallelCellLoop(GmCellMesh* m, bool activeOnly, T& taskObj,
370  int ntasks, CellPartitionStrategy strategy,
371  int nworkers, const QVector<int>& taskAffinity)
372 {
373  S_TRACE();
374 
375  // Add tasks to the manager.
376  bool ok;
377  if(!addCellTasks<GmCellMesh, T>(m, activeOnly, taskObj, ntasks, strategy, nworkers, taskAffinity, &ok))
378  {
379  if(!ok)
380  luaL_error(_threadManager->threadLuaEnv(0)->state(), "%s", luaPrintable(QObject::tr("Task execution aborted.")));
381  return;
382  }
383 
384  // Run tasks!
385  _threadManager->execTasks(nworkers);
386 }
387 
394 template <class M, class T>
395 bool GmTaskManager::addCellTasks(M* m, bool activeOnly, T& taskObj,
396  int ntasks, CellPartitionStrategy strategy, int nworkers,
397  const QVector<int>& taskAffinity, bool* ok)
398 {
399  S_TRACE();
400 
401  *ok = true;
402 
403  // Adjust execution parameters, converting possible default options to real values
404  adjustExecParameters(&ntasks, &nworkers);
405  if(strategy == CPS_DEFAULT)
406  strategy = _defCellStrategy;
407  assert(ntasks > 0 && nworkers >= 0 && strategy >= 0 && strategy < CPS_NUM_STRATEGIES);
408 
409  // Partition the mesh nodes
410  QVector<GmCellTaskDataSet*> taskData = buildCellPartition(m, ntasks, activeOnly, strategy);
411  if(taskData.isEmpty()) // No tasks
412  return false;
413 
414  // Check affinity vector
415  bool useAffinity = !taskAffinity.isEmpty();
416  if(useAffinity && taskData.size() != taskAffinity.size())
417  {
418  qDeleteAll(taskData);
419  gmErrorMsg(_threadManager->logger(), QObject::tr("Error creating cell tasks. Affinity vector size (%1) is different from the number of tasks (%2).")
420  .arg(taskAffinity.size()).arg(taskData.size()));
421  *ok = false;
422  return false; // No way to support the affinity request
423  }
424 
425  // Add tasks to the manager. The thread manager takes ownership and destroys the tasks
426  class CellTask: public GmThreadTask
427  {
428  public:
429  CellTask(T& task, GmCellTaskDataSet* dataSet) : _task(task), _dataSet(dataSet) {}
430  virtual ~CellTask() { delete _dataSet; }
431 
432  virtual GmThreadTaskResult run(const GmThread* thread) { return _task(_dataSet, thread); }
433 
434  private:
435  T& _task;
436  GmCellTaskDataSet* _dataSet;
437  };
438 
439  int i = 0;
440  foreach(GmCellTaskDataSet* ds, taskData)
441  _threadManager->addTask(new CellTask(taskObj, ds), useAffinity ? taskAffinity[i++] : 0);
442 
443  return true;
444 }
445 
446 #endif
void addTask(GmThreadTask *task, int tid=0)
Adds a task to the manager. The manager gets the task ownership and will delete it after execution....
Definition: gmThreadManager.cpp:400
The number of partition strategies.
Definition: gmTaskManager.h:107
bool addCellTasks(M *m, bool activeOnly, T &taskObj, int ntasks, CellPartitionStrategy strategy, int nworkers, const QVector< int > &taskAffinity, bool *ok)
Worker function for adding cell tasks to the thread manager. See runParallelCellLoop() for parameter ...
Definition: gmTaskManager.h:395
GmThreadManager * _threadManager
The thread manager used to run tasks in multiple threads.
Definition: gmTaskManager.h:163
Our thread wrapper, specialized for running tasks from the thread manager queue.
Definition: gmThreadManager.h:217
Base class for accessing the set of cells belonging to a mesh partition.
Definition: gmTaskManager.h:62
void execTasks(int nthreads)
Similar to runTasks() but raising a lua error if the execution was either canceled or aborted.
Definition: gmThreadManager.cpp:575
Interface for helping managing mesh cell groups and returning the elements in their union.
Definition: gmCellGroupSet.h:34
const GmLogCategory & logger()
Returns the thread manager logger.
Definition: gmThreadManager.h:147
GmAffectedNodes
Affected node types for node based values (GM_NODE_COORDINATES, GM_NODE_ATTRIBUTE or GM_NODE_STATEVAR...
Definition: gmValueInfo.h:91
Declaration of usefull configuration definitions for the Core library.
lua_State * state()
Declaration of the GmValueInfo class.
#define S_TRACE()
Macro for run time stack tracking at release build.
Definition: gmTrace.h:44
Base interface for mesh cells.
Definition: gmCell.h:81
GmThreadTaskResult
Possible results for running a task.
Definition: gmThreadManager.h:60
The thread finished its execution without errors.
Definition: gmThreadManager.h:62
void execParallelNodeLoop(const GmMesh *m, GmAffectedNodes affNodes, T &taskObj, int ntasks=0, NodePartitionStrategy strategy=NPS_DEFAULT, int nworkers=-1, const QVector< int > &taskAffinity=QVector< int >())
Similar to runParallelNodeLoop() but raising a lua error if the execution was either canceled or abor...
Definition: gmTaskManager.h:222
QVector< GmCellTaskDataSet * > buildCellPartition(M *m, int n, bool activeOnly, CellPartitionStrategy strategy) const
Creates a list with the 'n' parallel cell tasks to be executed, given the destination mesh,...
Definition: gmTaskManager.cpp:205
QString tr(const char *sourceText, const char *disambiguation, int n)
The number of partition strategies.
Definition: gmTaskManager.h:92
NodePartitionStrategy
Partition strategies for a parallel node loop.
Definition: gmTaskManager.h:81
bool addNodeTasks(const GmMesh *m, GmAffectedNodes affNodes, T &taskObj, int ntasks, NodePartitionStrategy strategy, int nworkers, const QVector< int > &taskAffinity, bool *ok)
Worker function for adding node tasks to the thread manager. See runParallelNodeLoop() for parameter ...
Definition: gmTaskManager.h:249
Auxiliar class used to store the complete set of simulation data.
Definition: gmSimulationData.h:51
Declaration of the GmThreadManager class.
int _defNumTasks
The default number of tasks for a parallel loop. Negative numbers are multiples of the number of work...
Definition: gmTaskManager.h:165
The thread aborted.
Definition: gmThreadManager.h:64
Base interface class for CellMesh type plugins.
Definition: gmCellMesh.h:39
Interface for a task executed by a thread manager thread.
Definition: gmThreadManager.h:68
GmThreadTaskResult runTasks(int nthreads)
Executes the set of tasks in the queue using at most 'nthreads'. Waits for the threads to finish befo...
Definition: gmThreadManager.cpp:444
LuaEnv * threadLuaEnv(int tid) const
Returns the Lua environment associated with the given thread id.
Definition: gmThreadManager.cpp:350
Uses the default strategy given by simulation options or CPS_SEQUENTIAL_BLOCK if missing.
Definition: gmTaskManager.h:98
Base class for accessing the set of nodes belonging to a mesh partition.
Definition: gmTaskManager.h:43
Uses the default strategy given by simulation options or NPS_SEQUENTIAL_BLOCK if missing.
Definition: gmTaskManager.h:83
virtual ~GmCellTaskDataSet()
Virtual destructor.
Definition: gmTaskManager.h:66
#define GMC_API_EXPORT
Macro for controling if the class is being exported (GEMA_CORE_LIB defined) or imported (GEMA_CORE_LI...
Definition: gmCoreConfig.h:35
GmThreadTaskResult runParallelNodeLoop(const GmMesh *m, GmAffectedNodes affNodes, T &taskObj, int ntasks=0, NodePartitionStrategy strategy=NPS_DEFAULT, int nworkers=-1, const QVector< int > &taskAffinity=QVector< int >())
Parallel executes the given task over mesh nodes.
Definition: gmTaskManager.h:205
bool isEmpty() const const
QVector< GmNodeTaskDataSet * > buildNodePartition(const GmMesh *m, int n, GmAffectedNodes affNodes, NodePartitionStrategy strategy) const
Creates a list with the 'n' parallel node tasks to be executed, given the destination mesh,...
Definition: gmTaskManager.cpp:159
int numTasks() const
Returns the number of tasks added to the queue by addNodetasks or addCellTasks.
Definition: gmTaskManager.h:119
virtual ~GmNodeTaskDataSet()
Virtual destructor.
Definition: gmTaskManager.h:47
CellPartitionStrategy _defCellStrategy
The default strategy for cell partitioning.
Definition: gmTaskManager.h:167
Thread manager used for handling parallel executions.
Definition: gmThreadManager.h:94
GmThreadManager * threadManager() const
Returns the associated thread manager.
Definition: gmTaskManager.h:116
NodePartitionStrategy _defNodeStrategy
The default strategy for node partitioning.
Definition: gmTaskManager.h:166
void adjustExecParameters(int *ntasks, int *nworkers)
Adjusts task execution parameters, converting possible default options to real values.
Definition: gmTaskManager.cpp:142
CellPartitionStrategy
Partition strategies for a parallel cell loop.
Definition: gmTaskManager.h:96
int size() const const
Base interface class for Mesh type plugins.
Definition: gmMesh.h:44
QString arg(qlonglong a, int fieldWidth, int base, QChar fillChar) const const
Auxiliary configuration file used to enable or disable compiling the GeMA tools with support for usin...
Task manager used for handling parallel executions.
Definition: gmTaskManager.h:76
void execParallelCellLoop(GmCellMesh *m, bool activeOnly, T &taskObj, int ntasks=0, CellPartitionStrategy strategy=CPS_DEFAULT, int nworkers=-1, const QVector< int > &taskAffinity=QVector< int >())
Similar to runParallelCellLoop() but raising a lua error if the execution was either canceled or abor...
Definition: gmTaskManager.h:369
GmThreadTaskResult runParallelCellLoop(GmCellMesh *m, bool activeOnly, T &taskObj, int ntasks=0, CellPartitionStrategy strategy=CPS_DEFAULT, int nworkers=-1, const QVector< int > &taskAffinity=QVector< int >())
Parallel executes the given task over mesh cells.
Definition: gmTaskManager.h:335