Arcane  v3.14.10.0
Documentation développeur
Chargement...
Recherche...
Aucune correspondance
TBBTaskImplementation.cc
1// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
2//-----------------------------------------------------------------------------
3// Copyright 2000-2023 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
4// See the top-level COPYRIGHT file for details.
5// SPDX-License-Identifier: Apache-2.0
6//-----------------------------------------------------------------------------
7/*---------------------------------------------------------------------------*/
8/* TBBTaskImplementation.cc (C) 2000-2023 */
9/* */
10/* Implémentation des tâches utilisant TBB (Intel Threads Building Blocks). */
11/*---------------------------------------------------------------------------*/
12/*---------------------------------------------------------------------------*/
13
14#include "arcane/utils/IThreadImplementation.h"
15#include "arcane/utils/NotImplementedException.h"
16#include "arcane/utils/IFunctor.h"
17#include "arcane/utils/CheckedConvert.h"
18#include "arcane/utils/ForLoopRanges.h"
20#include "arcane/utils/IObservable.h"
21#include "arcane/utils/PlatformUtils.h"
22#include "arcane/utils/Profiling.h"
23#include "arcane/utils/MemoryAllocator.h"
24
25#include "arcane/FactoryService.h"
26
27#include <new>
28#include <stack>
29
30// Il faut définir cette macro pour que la classe 'blocked_rangeNd' soit disponible
31
32#define TBB_PREVIEW_BLOCKED_RANGE_ND 1
33
34// la macro 'ARCANE_USE_ONETBB' est définie dans le CMakeLists.txt
35// si on compile avec la version OneTBB version 2021+
36// (https://github.com/oneapi-src/oneTBB.git)
37// A terme ce sera la seule version supportée par Arcane.
38
39#ifdef ARCANE_USE_ONETBB
40
41// Nécessaire pour avoir accès à task_scheduler_handle
42#define TBB_PREVIEW_WAITING_FOR_WORKERS 1
43#include <tbb/tbb.h>
44#include <oneapi/tbb/concurrent_set.h>
45#include <oneapi/tbb/global_control.h>
46
47#else // ARCANE_USE_ONETBB
48
49// NOTE GG: depuis mars 2019, la version 2018+ des TBB est obligatoire.
50#include <tbb/tbb.h>
51#include <tbb/blocked_rangeNd.h>
52/*
53 * Maintenant vérifie que la version est au moins 2018.
54 */
55#if (TBB_VERSION_MAJOR < 4) || (TBB_VERSION_MAJOR==4 && TBB_VERSION_MINOR<2)
56#define ARCANE_OLD_TBB
57#endif
58
59#if (TBB_VERSION_MAJOR<2018)
60#define ARCANE_OLD_TBB
61#endif
62
63#ifdef ARCANE_OLD_TBB
64#error "Your version of TBB is tool old. TBB 2018+ is required. Please disable TBB in configuration"
65#endif
66
67#include <thread>
68#include <mutex>
69
70#endif // ARCANE_USE_ONETBB
71
72/*---------------------------------------------------------------------------*/
73/*---------------------------------------------------------------------------*/
74
75namespace Arcane
76{
77
78class TBBTaskImplementation;
79
80// TODO: utiliser un pool mémoire spécifique pour gérer les
81// OneTBBTask pour optimiser les new/delete des instances de cette classe.
82// Auparavant avec les anciennes versions de TBB cela était géré avec
83// la méthode 'tbb::task::allocate_child()'.
84
85/*---------------------------------------------------------------------------*/
86/*---------------------------------------------------------------------------*/
87
88namespace
89{
90
91// Positif si on récupère les statistiques d'exécution
92bool isStatActive()
93{
95}
96
101class ScopedExecInfo
102{
103 public:
104
105 explicit ScopedExecInfo(const ForLoopRunInfo& run_info)
106 : m_run_info(run_info)
107 {
108 // Si run_info.execInfo() n'est pas nul, on l'utilise.
109 // Cela signifie que c'est l'appelant de qui va gérer les statistiques
110 // d'exécution. Sinon, on utilise \a m_stat_info si les statistiques
111 // d'exécution sont demandées.
112 ForLoopOneExecStat* ptr = run_info.execStat();
113 if (ptr){
114 m_stat_info_ptr = ptr;
115 m_use_own_run_info = false;
116 }
117 else
118 m_stat_info_ptr = isStatActive() ? &m_stat_info : nullptr;
119 }
120 ~ScopedExecInfo()
121 {
122 if (m_stat_info_ptr && m_use_own_run_info)
123 ProfilingRegistry::_threadLocalForLoopInstance()->merge(*m_stat_info_ptr,m_run_info.traceInfo());
124 }
125
126 public:
127
128 ForLoopOneExecStat* statInfo() const { return m_stat_info_ptr; }
129 bool isOwn() const { return m_use_own_run_info; }
130 private:
131
132 ForLoopOneExecStat m_stat_info;
133 ForLoopOneExecStat* m_stat_info_ptr = nullptr;
134 ForLoopRunInfo m_run_info;
136 bool m_use_own_run_info = true;
137};
138
139/*---------------------------------------------------------------------------*/
140/*---------------------------------------------------------------------------*/
141
142inline int _currentTaskTreadIndex()
143{
144 // NOTE: Avec OneTBB 2021, la valeur n'est plus '0' si on appelle cette méthode
145 // depuis un thread en dehors d'un task_arena. Avec la version 2021,
146 // la valeur est 65535.
147 // NOTE: Il semble que cela soit un bug de la 2021.3.
148 return tbb::this_task_arena::current_thread_index();
149}
150
151inline tbb::blocked_rangeNd<Int32,1>
152_toTBBRange(const ComplexForLoopRanges<1>& r)
153{
154 return {{r.lowerBound<0>(), r.upperBound<0>()}};
155}
156
157inline tbb::blocked_rangeNd<Int32,2>
158_toTBBRange(const ComplexForLoopRanges<2>& r)
159{
160 return {{r.lowerBound<0>(), r.upperBound<0>()},
161 {r.lowerBound<1>(), r.upperBound<1>()}};
162
163}
164
165inline tbb::blocked_rangeNd<Int32,3>
166_toTBBRange(const ComplexForLoopRanges<3>& r)
167{
168 return {{r.lowerBound<0>(), r.upperBound<0>()},
169 {r.lowerBound<1>(), r.upperBound<1>()},
170 {r.lowerBound<2>(), r.upperBound<2>()}};
171}
172
173inline tbb::blocked_rangeNd<Int32,4>
174_toTBBRange(const ComplexForLoopRanges<4>& r)
175{
176 return {{r.lowerBound<0>(), r.upperBound<0>()},
177 {r.lowerBound<1>(), r.upperBound<1>()},
178 {r.lowerBound<2>(), r.upperBound<2>()},
179 {r.lowerBound<3>(), r.upperBound<3>()}};
180}
181
182/*---------------------------------------------------------------------------*/
183/*---------------------------------------------------------------------------*/
184
185inline tbb::blocked_rangeNd<Int32,2>
186_toTBBRangeWithGrain(const tbb::blocked_rangeNd<Int32,2>& r,std::size_t grain_size)
187{
188 return {{r.dim(0).begin(), r.dim(0).end(), grain_size},
189 {r.dim(1).begin(), r.dim(1).end()}};
190}
191
192inline tbb::blocked_rangeNd<Int32,3>
193_toTBBRangeWithGrain(const tbb::blocked_rangeNd<Int32,3>& r,std::size_t grain_size)
194{
195 return {{r.dim(0).begin(), r.dim(0).end(), grain_size},
196 {r.dim(1).begin(), r.dim(0).end()},
197 {r.dim(2).begin(), r.dim(0).end()}};
198}
199
200inline tbb::blocked_rangeNd<Int32,4>
201_toTBBRangeWithGrain(const tbb::blocked_rangeNd<Int32,4>& r,std::size_t grain_size)
202{
203 return {{r.dim(0).begin(), r.dim(0).end(), grain_size},
204 {r.dim(1).begin(), r.dim(1).end()},
205 {r.dim(2).begin(), r.dim(2).end()},
206 {r.dim(3).begin(), r.dim(3).end()}};
207}
208
209/*---------------------------------------------------------------------------*/
210/*---------------------------------------------------------------------------*/
211
212inline ComplexForLoopRanges<2>
213_fromTBBRange(const tbb::blocked_rangeNd<Int32,2>& r)
214{
215 using BoundsType = ArrayBounds<MDDim2>;
216 using ArrayExtentType = typename BoundsType::ArrayExtentType;
217
218 BoundsType lower_bounds(ArrayExtentType(r.dim(0).begin(),r.dim(1).begin()));
219 auto s0 = static_cast<Int32>(r.dim(0).size());
220 auto s1 = static_cast<Int32>(r.dim(1).size());
221 BoundsType sizes(ArrayExtentType(s0,s1));
222 return { lower_bounds, sizes };
223}
224
225inline ComplexForLoopRanges<3>
226_fromTBBRange(const tbb::blocked_rangeNd<Int32,3> & r)
227{
228 using BoundsType = ArrayBounds<MDDim3>;
229 using ArrayExtentType = typename BoundsType::ArrayExtentType;
230
231 BoundsType lower_bounds(ArrayExtentType(r.dim(0).begin(),r.dim(1).begin(),r.dim(2).begin()));
232 auto s0 = static_cast<Int32>(r.dim(0).size());
233 auto s1 = static_cast<Int32>(r.dim(1).size());
234 auto s2 = static_cast<Int32>(r.dim(2).size());
235 BoundsType sizes(ArrayExtentType(s0,s1,s2));
236 return { lower_bounds, sizes };
237}
238
239inline ComplexForLoopRanges<4>
240_fromTBBRange(const tbb::blocked_rangeNd<Int32,4>& r)
241{
242 using BoundsType = ArrayBounds<MDDim4>;
243 using ArrayExtentType = typename BoundsType::ArrayExtentType;
244
245 BoundsType lower_bounds(ArrayExtentType(r.dim(0).begin(),r.dim(1).begin(),r.dim(2).begin(),r.dim(3).begin()));
246 auto s0 = static_cast<Int32>(r.dim(0).size());
247 auto s1 = static_cast<Int32>(r.dim(1).size());
248 auto s2 = static_cast<Int32>(r.dim(2).size());
249 auto s3 = static_cast<Int32>(r.dim(3).size());
250 BoundsType sizes(ArrayExtentType(s0,s1,s2,s3));
251 return { lower_bounds, sizes };
252}
253
254}
255
256/*---------------------------------------------------------------------------*/
257/*---------------------------------------------------------------------------*/
258
259#ifdef ARCANE_USE_ONETBB
260
261/*---------------------------------------------------------------------------*/
262/*---------------------------------------------------------------------------*/
263
264class OneTBBTaskFunctor
265{
266 public:
267 OneTBBTaskFunctor(ITaskFunctor* functor,ITask* task)
268 : m_functor(functor), m_task(task) {}
269 public:
270 void operator()() const
271 {
272 if (m_functor){
273 ITaskFunctor* tf = m_functor;
274 m_functor = 0;
275 TaskContext task_context(m_task);
276 //cerr << "FUNC=" << typeid(*tf).name();
277 tf->executeFunctor(task_context);
278 }
279 }
280 public:
281 mutable ITaskFunctor* m_functor;
282 ITask* m_task;
283};
284
285/*---------------------------------------------------------------------------*/
286/*---------------------------------------------------------------------------*/
287
288class OneTBBTask
289: public ITask
290{
291 public:
292 static const int FUNCTOR_CLASS_SIZE = 32;
293 public:
294 OneTBBTask(ITaskFunctor* f)
295 : m_functor(f)
296 {
297 m_functor = f->clone(functor_buf,FUNCTOR_CLASS_SIZE);
298 }
299 public:
300 OneTBBTaskFunctor taskFunctor() { return OneTBBTaskFunctor(m_functor,this); }
301 void launchAndWait() override;
302 void launchAndWait(ConstArrayView<ITask*> tasks) override;
303 protected:
304 virtual ITask* _createChildTask(ITaskFunctor* functor) override;
305 public:
306 ITaskFunctor* m_functor;
307 char functor_buf[FUNCTOR_CLASS_SIZE];
308};
309using TBBTask = OneTBBTask;
310
311/*---------------------------------------------------------------------------*/
312/*---------------------------------------------------------------------------*/
313
314#else // ARCANE_USE_ONETBB
315
316/*---------------------------------------------------------------------------*/
317/*---------------------------------------------------------------------------*/
318
320: public tbb::task
321, public ITask
322{
323 public:
324 static const int FUNCTOR_CLASS_SIZE = 32;
325 public:
327 : m_functor(f)
328 {
329 m_functor = f->clone(functor_buf,FUNCTOR_CLASS_SIZE);
330 }
331 public:
332 tbb::task* execute() override
333 {
334 if (m_functor){
335 ITaskFunctor* tf = m_functor;
336 m_functor = 0;
338 //cerr << "FUNC=" << typeid(*tf).name();
339 tf->executeFunctor(task_context);
340 }
341 return nullptr;
342 }
343 void launchAndWait() override;
345 protected:
346 ITask* _createChildTask(ITaskFunctor* functor) final;
347 public:
348 ITaskFunctor* m_functor;
349 char functor_buf[FUNCTOR_CLASS_SIZE];
350};
351using TBBTask = LegacyTBBTask;
352
353/*---------------------------------------------------------------------------*/
354/*---------------------------------------------------------------------------*/
355
356#endif // ARCANE_USE_ONETBB
357
358/*---------------------------------------------------------------------------*/
359/*---------------------------------------------------------------------------*/
360/*
361 * Ne pas utiliser l'observer locale au task_arena.
362 * Utiliser l'observer global au scheduler.
363 * Pour l'id, utiliser tbb::this_task_arena::current_thread_index().
364 */
366: public ITaskImplementation
367{
368 class Impl;
369 class ParallelForExecute;
370 template<int RankValue>
372
373 public:
374 // Pour des raisons de performance, s'aligne sur une ligne de cache
375 // et utilise un padding.
377 {
378 public:
379 TaskThreadInfo() : m_task_index(-1){}
380 public:
381 void setTaskIndex(Integer v) { m_task_index = v; }
382 Integer taskIndex() const { return m_task_index; }
383 private:
384 Integer m_task_index;
385 };
394 {
395 public:
397 : m_tti(tti), m_old_task_index(-1)
398 {
399 if (tti){
400 m_old_task_index = tti->taskIndex();
401 tti->setTaskIndex(task_index);
402 }
403 }
405 {
406 if (m_tti)
407 m_tti->setTaskIndex(m_old_task_index);
408 }
409 private:
410 TaskThreadInfo* m_tti;
411 Integer m_old_task_index;
412 };
413
414 public:
416 : m_is_active(false), m_p(nullptr)
417 {
418 ARCANE_UNUSED(sbi);
419 }
420 ~TBBTaskImplementation() override;
421 public:
422 void build() {}
423 void initialize(Int32 nb_thread) override;
424 void terminate() override;
425
427 {
428#ifdef ARCANE_USE_ONETBB
429 OneTBBTask* t = new OneTBBTask(f);
430#else
431 LegacyTBBTask* t = new(tbb::task::allocate_root()) LegacyTBBTask(f);
432#endif
433 return t;
434 }
435
436 void executeParallelFor(Int32 begin,Int32 size,const ParallelLoopOptions& options,IRangeFunctor* f) final;
437 void executeParallelFor(Int32 begin,Int32 size,Integer grain_size,IRangeFunctor* f) final;
438 void executeParallelFor(Int32 begin,Int32 size,IRangeFunctor* f) final
439 {
440 executeParallelFor(begin,size,TaskFactory::defaultParallelLoopOptions(),f);
441 }
442 void executeParallelFor(const ParallelFor1DLoopInfo& loop_info) override;
443
445 const ParallelLoopOptions& options,
446 IMDRangeFunctor<1>* functor) final
447 {
448 _executeMDParallelFor<1>(loop_ranges,functor,options);
449 }
451 const ParallelLoopOptions& options,
452 IMDRangeFunctor<2>* functor) final
453 {
454 _executeMDParallelFor<2>(loop_ranges,functor,options);
455 }
457 const ParallelLoopOptions& options,
458 IMDRangeFunctor<3>* functor) final
459 {
460 _executeMDParallelFor<3>(loop_ranges,functor,options);
461 }
463 const ParallelLoopOptions& options,
464 IMDRangeFunctor<4>* functor) final
465 {
466 _executeMDParallelFor<4>(loop_ranges,functor,options);
467 }
468
469 bool isActive() const final
470 {
471 return m_is_active;
472 }
473
474 Int32 nbAllowedThread() const final;
475
477 {
478 return (nbAllowedThread() <= 1 ) ? 0 : _currentTaskTreadIndex();
479 }
480
481 Int32 currentTaskIndex() const final;
482
483 void printInfos(std::ostream& o) const final;
484
485 public:
486
493 TaskThreadInfo* currentTaskThreadInfo() const;
494
495 private:
496
497 bool m_is_active;
498 Impl* m_p;
499
500 private:
501
502 template<int RankValue> void
504 IMDRangeFunctor<RankValue>* functor,
505 const ParallelLoopOptions& options);
506 void _executeParallelFor(const ParallelFor1DLoopInfo& loop_info);
507};
508
509/*---------------------------------------------------------------------------*/
510/*---------------------------------------------------------------------------*/
511
513{
515 : public tbb::task_scheduler_observer
516 {
517 public:
519 :
520#ifdef ARCANE_USE_ONETBB
521 tbb::task_scheduler_observer(p->m_main_arena),
522#endif
523 m_p(p)
524 {
525 }
526 void on_scheduler_entry(bool is_worker) override
527 {
528 m_p->notifyThreadCreated(is_worker);
529 }
530 void on_scheduler_exit(bool is_worker) override
531 {
532 m_p->notifyThreadDestroyed(is_worker);
533 }
535 };
536
537 public:
538 Impl() :
539 m_task_observer(this),
540 m_thread_task_infos(AlignedMemoryAllocator::CacheLine())
541 {
542#ifdef ARCANE_USE_ONETBB
543 m_nb_allowed_thread = tbb::info::default_concurrency();
544#else
545 m_nb_allowed_thread = tbb::task_scheduler_init::default_num_threads();
546#endif
547 _init();
548 }
549 Impl(Int32 nb_thread)
550 :
552 m_scheduler_init(nb_thread),
553#endif
554 m_main_arena(nb_thread),
555 m_task_observer(this),
556 m_thread_task_infos(AlignedMemoryAllocator::CacheLine())
557 {
558 m_nb_allowed_thread = nb_thread;
559 _init();
560 }
561 public:
562 Int32 nbAllowedThread() const { return m_nb_allowed_thread; }
563 TaskThreadInfo* threadTaskInfo(Integer index) { return &m_thread_task_infos[index]; }
564 private:
565 Int32 m_nb_allowed_thread;
566
567 public:
568 void terminate()
569 {
570 for( auto x : m_sub_arena_list ){
571 if (x)
572 x->terminate();
573 delete x;
574 }
575 m_sub_arena_list.clear();
576 m_main_arena.terminate();
577#ifdef ARCANE_USE_ONETBB
578 m_task_observer.observe(false);
579 oneapi::tbb::finalize(m_task_scheduler_handle);
580#else
581 m_scheduler_init.terminate();
582 m_task_observer.observe(false);
583#endif
584 }
585 public:
586 void notifyThreadCreated(bool is_worker)
587 {
588 std::thread::id my_thread_id = std::this_thread::get_id();
589
590#ifdef ARCANE_USE_ONETBB
591 // Avec OneTBB, cette méthode est appelée à chaque fois qu'on rentre
592 // dans notre 'task_arena'. Comme il ne faut appeler qu'une seule
593 // fois la méthode de notification on utilise un ensemble pour
594 // conserver la liste des threads déjà créés.
595 // NOTE: On ne peut pas utiliser cette méthode avec la version TBB historique
596 // (2018) car cette méthode 'contains' n'existe pas
597 if (m_constructed_thread_map.contains(my_thread_id))
598 return;
599 m_constructed_thread_map.insert(my_thread_id);
600#endif
601
602 // Il faut toujours un verrou car on n'est pas certain que
603 // les méthodes appelées par l'observable soient thread-safe
604 // (et aussi TaskFactory::createThreadObservable() ne l'est pas)
605 {
606 std::scoped_lock sl(m_thread_created_mutex);
608 std::cout << "TBB: CREATE THREAD"
609 << " nb_allowed=" << m_nb_allowed_thread
610#ifdef ARCANE_USE_ONETBB
611 << " tbb_default_allowed=" << tbb::info::default_concurrency()
612#else
613 << " tbb_default_allowed=" << tbb::task_scheduler_init::default_num_threads()
614#endif
615 << " id=" << my_thread_id
616 << " arena_id=" << _currentTaskTreadIndex()
617 << " is_worker=" << is_worker
618 << "\n";
619 }
620 TaskFactory::createThreadObservable()->notifyAllObservers();
621 }
622 }
623
624 void notifyThreadDestroyed([[maybe_unused]] bool is_worker)
625 {
626#ifdef ARCANE_USE_ONETBB
627 // Avec OneTBB, cette méthode est appelée à chaque fois qu'on sort
628 // de l'arène principale. Du coup elle ne correspond pas vraiment à une
629 // destruction de thread. On ne fait donc rien pour cette notification
630 // TODO: Regarder comment on peut être notifié de la destruction effective
631 // du thread.
632#else
633 // Il faut toujours un verrou car on n'est pas certain que
634 // les méthodes appelées par l'observable soient thread-safe
635 // (et aussi TaskFactory::createThreadObservable() ne l'est pas)
636 std::scoped_lock sl(m_thread_created_mutex);
638 std::cout << "TBB: DESTROY THREAD"
639 << " id=" << std::this_thread::get_id()
640 << " arena_id=" << _currentTaskTreadIndex()
641 << " is_worker=" << is_worker
642 << '\n';
643 }
644 TaskFactory::destroyThreadObservable()->notifyAllObservers();
645#endif
646 }
647 private:
648#ifdef ARCANE_USE_ONETBB
649#if TBB_VERSION_MAJOR>2021 || (TBB_VERSION_MAJOR==2021 && TBB_VERSION_MINOR>5)
650 oneapi::tbb::task_scheduler_handle m_task_scheduler_handle = oneapi::tbb::attach();
651#else
652 oneapi::tbb::task_scheduler_handle m_task_scheduler_handle = tbb::task_scheduler_handle::get();
653#endif
654#else
655 tbb::task_scheduler_init m_scheduler_init;
656#endif
657 public:
658 tbb::task_arena m_main_arena;
660 std::vector<tbb::task_arena*> m_sub_arena_list;
661 private:
662 TaskObserver m_task_observer;
663 std::mutex m_thread_created_mutex;
664 UniqueArray<TaskThreadInfo> m_thread_task_infos;
665#ifdef ARCANE_USE_ONETBB
666 tbb::concurrent_set<std::thread::id> m_constructed_thread_map;
667#endif
668 void _init()
669 {
671 std::cout << "TBB: TBBTaskImplementationInit nb_allowed_thread=" << m_nb_allowed_thread
672 << " id=" << std::this_thread::get_id()
673 << " version=" << TBB_VERSION_MAJOR << "." << TBB_VERSION_MINOR
674 << "\n";
675 }
676 m_thread_task_infos.resize(m_nb_allowed_thread);
677 m_task_observer.observe(true);
678 Integer max_arena_size = m_nb_allowed_thread;
679 // Limite artificiellement le nombre de tbb::task_arena
680 // pour éviter d'avoir trop d'objets alloués.
681 if (max_arena_size>512)
682 max_arena_size = 512;
683 m_sub_arena_list.resize(max_arena_size);
684 m_sub_arena_list[0] = m_sub_arena_list[1] = nullptr;
685 for( Integer i=2; i<max_arena_size; ++i )
686 m_sub_arena_list[i] = new tbb::task_arena(i);
687 }
688};
689
690/*---------------------------------------------------------------------------*/
691/*---------------------------------------------------------------------------*/
696{
697 public:
699 : m_functor(f), m_stat_info(stat_info), m_nb_allowed_thread(nb_allowed_thread){}
700 public:
701
702 void operator()(tbb::blocked_range<Integer>& range) const
703 {
704#ifdef ARCANE_CHECK
705 if (TaskFactory::verboseLevel()>=3){
706 std::ostringstream o;
707 o << "TBB: INDEX=" << TaskFactory::currentTaskThreadIndex()
708 << " id=" << std::this_thread::get_id()
709 << " max_allowed=" << m_nb_allowed_thread
710 << " range_begin=" << range.begin() << " range_size=" << range.size()
711 << "\n";
712 std::cout << o.str();
713 std::cout.flush();
714 }
715
717 if (tbb_index<0 || tbb_index>=m_nb_allowed_thread)
718 ARCANE_FATAL("Invalid index for thread idx={0} valid_interval=[0..{1}[",
719 tbb_index,m_nb_allowed_thread);
720#endif
721
722 if (m_stat_info)
723 m_stat_info->incrementNbChunk();
724 m_functor->executeFunctor(range.begin(),CheckedConvert::toInteger(range.size()));
725 }
726
727 private:
728 IRangeFunctor* m_functor;
729 ForLoopOneExecStat* m_stat_info = nullptr;
730 Int32 m_nb_allowed_thread;
731};
732
733/*---------------------------------------------------------------------------*/
734/*---------------------------------------------------------------------------*/
738template<int RankValue>
740{
741 public:
742
744 : m_functor(f), m_stat_info(stat_info), m_nb_allowed_thread(nb_allowed_thread){}
745
746 public:
747
748 void operator()(tbb::blocked_rangeNd<Int32,RankValue>& range) const
749 {
750#ifdef ARCANE_CHECK
751 if (TaskFactory::verboseLevel()>=3){
752 std::ostringstream o;
753 o << "TBB: INDEX=" << TaskFactory::currentTaskThreadIndex()
754 << " id=" << std::this_thread::get_id()
755 << " max_allowed=" << m_nb_allowed_thread
756 //<< " range_begin=" << range.begin() << " range_size=" << range.size()
757 << "\n";
758 std::cout << o.str();
759 std::cout.flush();
760 }
761
763 if (tbb_index<0 || tbb_index>=m_nb_allowed_thread)
764 ARCANE_FATAL("Invalid index for thread idx={0} valid_interval=[0..{1}[",
765 tbb_index,m_nb_allowed_thread);
766#endif
767
768 if (m_stat_info)
769 m_stat_info->incrementNbChunk();
770 m_functor->executeFunctor(_fromTBBRange(range));
771 }
772
773 private:
774
775 IMDRangeFunctor<RankValue>* m_functor = nullptr;
776 ForLoopOneExecStat* m_stat_info = nullptr;
777 Int32 m_nb_allowed_thread;
778};
779
780/*---------------------------------------------------------------------------*/
781/*---------------------------------------------------------------------------*/
800{
801 public:
803 Integer begin_index,Integer size,Integer grain_size,Integer nb_thread)
804 : m_impl(impl), m_tbb_for(tbb_for), m_nb_thread(nb_thread), m_begin_index(begin_index), m_size(size),
805 m_grain_size(grain_size), m_nb_block(0), m_block_size(0), m_nb_block_per_thread(0)
806 {
807 if (m_nb_thread<1)
808 m_nb_thread = 1;
809
810 if (m_grain_size>0){
811 m_block_size = m_grain_size;
812 if (m_block_size>0){
813 m_nb_block = m_size / m_block_size;
814 if ((m_size % m_block_size)!=0)
815 ++m_nb_block;
816 }
817 else
818 m_nb_block = 1;
819 m_nb_block_per_thread = m_nb_block / m_nb_thread;
820 if ((m_nb_block % m_nb_thread) != 0)
821 ++m_nb_block_per_thread;
822 }
823 else{
824 if (m_nb_block<1)
825 m_nb_block = m_nb_thread;
826 m_block_size = m_size / m_nb_block;
827 m_nb_block_per_thread = 1;
828 }
829 if (TaskFactory::verboseLevel()>=2){
830 std::cout << "TBBDeterministicParallelFor: BEGIN=" << m_begin_index << " size=" << m_size
831 << " grain_size=" << m_grain_size
832 << " nb_block=" << m_nb_block << " nb_thread=" << m_nb_thread
833 << " nb_block_per_thread=" << m_nb_block_per_thread
834 << " block_size=" << m_block_size
835 << " block_size*nb_block=" << m_block_size*m_nb_block << '\n';
836 }
837 }
838 public:
839
846 void operator()(tbb::blocked_range<Integer>& range) const
847 {
848 Integer nb_iter = static_cast<Integer>(range.size());
849 for( Integer i=0; i<nb_iter; ++i ){
850 Integer task_id = range.begin() + i;
851 for ( Integer k=0, kn=m_nb_block_per_thread; k<kn; ++k ){
852 Integer block_id = task_id + (k * m_nb_thread);
853 if (block_id<m_nb_block)
854 _doBlock(task_id,block_id);
855 }
856 }
857 }
858
859 void _doBlock(Integer task_id,Integer block_id) const
860 {
861 TBBTaskImplementation::TaskInfoLockGuard guard(m_impl->currentTaskThreadInfo(),task_id);
862
863 Integer iter_begin = block_id * m_block_size;
864 Integer iter_size = m_block_size;
865 if ((block_id+1)==m_nb_block){
866 // Pour le dernier bloc, la taille est le nombre d'éléments restants
867 iter_size = m_size - iter_begin;
868 }
869 iter_begin += m_begin_index;
870#ifdef ARCANE_CHECK
871 if (TaskFactory::verboseLevel()>=3){
872 std::ostringstream o;
873 o << "TBB: DoBlock: BLOCK task_id=" << task_id << " block_id=" << block_id
874 << " iter_begin=" << iter_begin << " iter_size=" << iter_size << '\n';
875 std::cout << o.str();
876 std::cout.flush();
877 }
878#endif
879 if (iter_size>0){
880 auto r = tbb::blocked_range<int>(iter_begin,iter_begin + iter_size);
881 m_tbb_for(r);
882 }
883 }
884
885 private:
886
887 TBBTaskImplementation* m_impl;
888 const TBBParallelFor& m_tbb_for;
889 Integer m_nb_thread;
890 Integer m_begin_index;
891 Integer m_size;
892 Integer m_grain_size;
893 Integer m_nb_block;
894 Integer m_block_size;
895 Integer m_nb_block_per_thread;
896};
897
898/*---------------------------------------------------------------------------*/
899/*---------------------------------------------------------------------------*/
900
902{
903 public:
904
906 Integer begin,Integer size,IRangeFunctor* f,ForLoopOneExecStat* stat_info)
907 : m_impl(impl), m_begin(begin), m_size(size), m_functor(f), m_options(options), m_stat_info(stat_info){}
908
909 public:
910
911 void operator()() const
912 {
913 Integer nb_thread = m_options.maxThread();
914 TBBParallelFor pf(m_functor,nb_thread,m_stat_info);
915 Integer gsize = m_options.grainSize();
916 tbb::blocked_range<Integer> range(m_begin,m_begin+m_size);
917 if (TaskFactory::verboseLevel()>=1)
918 std::cout << "TBB: TBBTaskImplementationInit ParallelForExecute begin=" << m_begin
919 << " size=" << m_size << " gsize=" << gsize
920 << " partitioner=" << (int)m_options.partitioner()
921 << " nb_thread=" << nb_thread << '\n';
922
923 if (gsize>0)
924 range = tbb::blocked_range<Integer>(m_begin,m_begin+m_size,gsize);
925
926 if (m_options.partitioner()==ParallelLoopOptions::Partitioner::Static){
927 tbb::parallel_for(range,pf,tbb::static_partitioner());
928 }
929 else if (m_options.partitioner()==ParallelLoopOptions::Partitioner::Deterministic){
930 tbb::blocked_range<Integer> range2(0,nb_thread,1);
931 TBBDeterministicParallelFor dpf(m_impl,pf,m_begin,m_size,gsize,nb_thread);
932 tbb::parallel_for(range2,dpf);
933 }
934 else
935 tbb::parallel_for(range,pf);
936 }
937 private:
938 TBBTaskImplementation* m_impl = nullptr;
939 Integer m_begin;
940 Integer m_size;
941 IRangeFunctor* m_functor = nullptr;
942 ParallelLoopOptions m_options;
943 ForLoopOneExecStat* m_stat_info = nullptr;
944};
945
946/*---------------------------------------------------------------------------*/
947/*---------------------------------------------------------------------------*/
948
949template<int RankValue>
951{
952 public:
953
955 const ParallelLoopOptions& options,
958 : m_impl(impl), m_tbb_range(_toTBBRange(range)), m_functor(f), m_options(options)
959 {
960 // On ne peut pas modifier les valeurs d'une instance de tbb::blocked_rangeNd.
961 // Il faut donc en reconstruire une complètement.
962
963 Integer gsize = m_options.grainSize();
964 if (gsize>0){
965 // Modifie la taille du grain pour la première dimension.
966 // TODO: pouvoir aussi modifier la valeur de 'grain_size' pour les autres dimensions.
967 m_tbb_range = _toTBBRangeWithGrain(m_tbb_range,gsize);
968 }
969 }
970
971 public:
972
973 void operator()() const
974 {
975 Integer nb_thread = m_options.maxThread();
976 TBBMDParallelFor<RankValue> pf(m_functor,nb_thread,m_stat_info);
977
978 if (m_options.partitioner()==ParallelLoopOptions::Partitioner::Static){
979 tbb::parallel_for(m_tbb_range,pf,tbb::static_partitioner());
980 }
981 else if (m_options.partitioner()==ParallelLoopOptions::Partitioner::Deterministic){
982 // TODO: implémenter le mode déterministe
983 ARCANE_THROW(NotImplementedException,"ParallelLoopOptions::Partitioner::Deterministic for multi-dimensionnal loops");
984 //tbb::blocked_range<Integer> range2(0,nb_thread,1);
985 //TBBDeterministicParallelFor dpf(m_impl,pf,m_begin,m_size,gsize,nb_thread);
986 //tbb::parallel_for(range2,dpf);
987 }
988 else
989 tbb::parallel_for(m_tbb_range,pf);
990 }
991 private:
992 TBBTaskImplementation* m_impl = nullptr;
993 tbb::blocked_rangeNd<Int32,RankValue> m_tbb_range;
994 IMDRangeFunctor<RankValue>* m_functor = nullptr;
995 ParallelLoopOptions m_options;
996 ForLoopOneExecStat* m_stat_info = nullptr;
997};
998
999/*---------------------------------------------------------------------------*/
1000/*---------------------------------------------------------------------------*/
1001
1002TBBTaskImplementation::
1003~TBBTaskImplementation()
1004{
1005 delete m_p;
1006}
1007
1008/*---------------------------------------------------------------------------*/
1009/*---------------------------------------------------------------------------*/
1010
1012initialize(Int32 nb_thread)
1013{
1014 if (nb_thread<0)
1015 nb_thread = 0;
1016 m_is_active = (nb_thread!=1);
1017 if (nb_thread!=0)
1018 m_p = new Impl(nb_thread);
1019 else
1020 m_p = new Impl();
1022 opts.setMaxThread(nbAllowedThread());
1024}
1025
1026/*---------------------------------------------------------------------------*/
1027/*---------------------------------------------------------------------------*/
1028
1030terminate()
1031{
1032 m_p->terminate();
1033}
1034
1035/*---------------------------------------------------------------------------*/
1036/*---------------------------------------------------------------------------*/
1037
1039nbAllowedThread() const
1040{
1041 return m_p->nbAllowedThread();
1042}
1043
1044/*---------------------------------------------------------------------------*/
1045/*---------------------------------------------------------------------------*/
1046
1048printInfos(std::ostream& o) const
1049{
1050#ifdef ARCANE_USE_ONETBB
1051 o << "OneTBBTaskImplementation"
1052 << " version=" << TBB_VERSION_STRING
1053 << " interface=" << TBB_INTERFACE_VERSION
1054 << " runtime_interface=" << TBB_runtime_interface_version();
1055#else
1056 o << "TBBTaskImplementation"
1057 << " version=" << TBB_VERSION_MAJOR << "." << TBB_VERSION_MINOR
1058 << " interface=" << TBB_INTERFACE_VERSION;
1059#endif
1060}
1061
1062/*---------------------------------------------------------------------------*/
1063/*---------------------------------------------------------------------------*/
1064
1065void TBBTaskImplementation::
1066_executeParallelFor(const ParallelFor1DLoopInfo& loop_info)
1067{
1068 ScopedExecInfo sei(loop_info.runInfo());
1069 ForLoopOneExecStat* stat_info = sei.statInfo();
1071
1072 Int32 begin = loop_info.beginIndex();
1073 Int32 size = loop_info.size();
1074 ParallelLoopOptions options = loop_info.runInfo().options().value_or(TaskFactory::defaultParallelLoopOptions());
1075 IRangeFunctor* f = loop_info.functor();
1076
1077 Integer max_thread = options.maxThread();
1078 Integer nb_allowed_thread = m_p->nbAllowedThread();
1079 if (max_thread<0)
1081
1083 std::cout << "TBB: TBBTaskImplementation executeParallelFor begin=" << begin
1084 << " size=" << size << " max_thread=" << max_thread
1085 << " grain_size=" << options.grainSize()
1086 << " nb_allowed=" << nb_allowed_thread << '\n';
1087
1088 // En exécution séquentielle, appelle directement la méthode \a f.
1089 if (max_thread==1 || max_thread==0){
1090 f->executeFunctor(begin,size);
1091 return;
1092 }
1093
1094 // Remplace les valeurs non initialisées de \a options par celles de \a m_default_loop_options
1095 ParallelLoopOptions true_options(options);
1096 true_options.mergeUnsetValues(TaskFactory::defaultParallelLoopOptions());
1097 true_options.setMaxThread(max_thread);
1098
1099 ParallelForExecute pfe(this,true_options,begin,size,f,stat_info);
1100
1101 tbb::task_arena* used_arena = nullptr;
1102 if (max_thread<nb_allowed_thread)
1103 used_arena = m_p->m_sub_arena_list[max_thread];
1104 if (!used_arena)
1105 used_arena = &(m_p->m_main_arena);
1106 used_arena->execute(pfe);
1107}
1108
1109/*---------------------------------------------------------------------------*/
1110/*---------------------------------------------------------------------------*/
1111
1112void TBBTaskImplementation::
1113executeParallelFor(const ParallelFor1DLoopInfo& loop_info)
1114{
1115 _executeParallelFor(loop_info);
1116}
1117
1118/*---------------------------------------------------------------------------*/
1119/*---------------------------------------------------------------------------*/
1125template<int RankValue> void TBBTaskImplementation::
1128 const ParallelLoopOptions& options)
1129{
1130 ScopedExecInfo sei(ForLoopRunInfo{});
1131 ForLoopOneExecStat* stat_info = sei.statInfo();
1132 impl::ScopedStatLoop scoped_loop(sei.isOwn() ? stat_info : nullptr);
1133
1135 std::cout << "TBB: TBBTaskImplementation executeMDParallelFor nb_dim=" << RankValue << '\n';
1136 Integer max_thread = options.maxThread();
1137 // En exécution séquentielle, appelle directement la méthode \a f.
1138 if (max_thread==1 || max_thread==0){
1139 functor->executeFunctor(loop_ranges);
1140 return;
1141 }
1142
1143 // Remplace les valeurs non initialisées de \a options par celles de \a m_default_loop_options
1146
1147 Integer nb_allowed_thread = m_p->nbAllowedThread();
1148 if (max_thread<0)
1150 tbb::task_arena* used_arena = nullptr;
1152 used_arena = m_p->m_sub_arena_list[max_thread];
1153 if (!used_arena)
1154 used_arena = &(m_p->m_main_arena);
1155
1156 // Pour l'instant pour la dimension 1, utilise le 'ParallelForExecute' historique
1157 if constexpr (RankValue==1){
1159 auto x1 = [&](Integer begin,Integer size)
1160 {
1161 functor->executeFunctor(makeLoopRanges(ForLoopRange(begin,size)));
1162 //functor->executeFunctor(ComplexForLoopRanges<1>(begin,size));
1163 };
1164 LambdaRangeFunctorT<decltype(x1)> functor_1d(x1);
1165 Integer begin1 = CheckedConvert::toInteger(range_1d.dim(0).begin());
1166 Integer size1 = CheckedConvert::toInteger(range_1d.dim(0).size());
1168 used_arena->execute(pfe);
1169 }
1170 else{
1172 used_arena->execute(pfe);
1173 }
1174}
1175
1176/*---------------------------------------------------------------------------*/
1177/*---------------------------------------------------------------------------*/
1178
1179void TBBTaskImplementation::
1180executeParallelFor(Integer begin,Integer size,Integer grain_size,IRangeFunctor* f)
1181{
1183 opts.setGrainSize(grain_size);
1185 executeParallelFor(ParallelFor1DLoopInfo(begin,size,f,run_info));
1186}
1187
1188/*---------------------------------------------------------------------------*/
1189/*---------------------------------------------------------------------------*/
1190
1191void TBBTaskImplementation::
1192executeParallelFor(Integer begin,Integer size,const ParallelLoopOptions& options,IRangeFunctor* f)
1193{
1194 executeParallelFor(ParallelFor1DLoopInfo(begin,size,f,ForLoopRunInfo(options)));
1195}
1196
1197/*---------------------------------------------------------------------------*/
1198/*---------------------------------------------------------------------------*/
1199
1202{
1204 if (thread_id>=0)
1205 return m_p->threadTaskInfo(thread_id);
1206 return nullptr;
1207}
1208
1209/*---------------------------------------------------------------------------*/
1210/*---------------------------------------------------------------------------*/
1211
1213currentTaskIndex() const
1214{
1216#ifdef ARCANE_USE_ONETBB
1217 if (thread_id<0 || thread_id>=m_p->nbAllowedThread())
1218 return 0;
1219#endif
1221 if (tti){
1222 Int32 task_index = tti->taskIndex();
1223 if (task_index>=0)
1224 return task_index;
1225 }
1226 return thread_id;
1227}
1228
1229/*---------------------------------------------------------------------------*/
1230/*---------------------------------------------------------------------------*/
1231
1232#ifdef ARCANE_USE_ONETBB
1233
1234/*---------------------------------------------------------------------------*/
1235/*---------------------------------------------------------------------------*/
1236
1237void OneTBBTask::
1238launchAndWait()
1239{
1240 tbb::task_group task_group;
1241 task_group.run(taskFunctor());
1242 task_group.wait();
1243 delete this;
1244}
1245
1246/*---------------------------------------------------------------------------*/
1247/*---------------------------------------------------------------------------*/
1248
1249void OneTBBTask::
1250launchAndWait(ConstArrayView<ITask*> tasks)
1251{
1252 tbb::task_group task_group;
1253 Integer n = tasks.size();
1254 if (n==0)
1255 return;
1256
1257 //set_ref_count(n+1);
1258 for( Integer i=0; i<n; ++i ){
1259 OneTBBTask* t = static_cast<OneTBBTask*>(tasks[i]);
1260 task_group.run(t->taskFunctor());
1261 }
1262 task_group.wait();
1263 for( Integer i=0; i<n; ++i ){
1264 OneTBBTask* t = static_cast<OneTBBTask*>(tasks[i]);
1265 delete t;
1266 }
1267}
1268
1269/*---------------------------------------------------------------------------*/
1270/*---------------------------------------------------------------------------*/
1271
1272ITask* OneTBBTask::
1273_createChildTask(ITaskFunctor* functor)
1274{
1275 OneTBBTask* t = new OneTBBTask(functor);
1276 return t;
1277}
1278
1279/*---------------------------------------------------------------------------*/
1280/*---------------------------------------------------------------------------*/
1281
1282#else // ARCANE_USE_ONETBB
1283
1284/*---------------------------------------------------------------------------*/
1285/*---------------------------------------------------------------------------*/
1286
1289{
1290 task::spawn_root_and_wait(*this);
1291}
1292
1293/*---------------------------------------------------------------------------*/
1294/*---------------------------------------------------------------------------*/
1295
1298{
1299 Integer n = tasks.size();
1300 if (n==0)
1301 return;
1302
1303 set_ref_count(n+1);
1304 for( Integer i=0; i<n-1; ++i ){
1305 TBBTask* t = static_cast<TBBTask*>(tasks[i]);
1306 spawn(*t);
1307 }
1308 spawn_and_wait_for_all(*static_cast<TBBTask*>(tasks[n-1]));
1309}
1310
1311/*---------------------------------------------------------------------------*/
1312/*---------------------------------------------------------------------------*/
1313
1314ITask* LegacyTBBTask::
1315_createChildTask(ITaskFunctor* functor)
1316{
1317 TBBTask* t = new(allocate_child()) TBBTask(functor);
1318 return t;
1319}
1320
1321/*---------------------------------------------------------------------------*/
1322/*---------------------------------------------------------------------------*/
1323
1324#endif // ARCANE_USE_ONETBB
1325
1326/*---------------------------------------------------------------------------*/
1327/*---------------------------------------------------------------------------*/
1328
1329ARCANE_REGISTER_APPLICATION_FACTORY(TBBTaskImplementation,ITaskImplementation,
1330 TBBTaskImplementation);
1331
1332/*---------------------------------------------------------------------------*/
1333/*---------------------------------------------------------------------------*/
1334
1335} // End namespace Arcane
1336
1337/*---------------------------------------------------------------------------*/
1338/*---------------------------------------------------------------------------*/
#define ARCANE_THROW(exception_class,...)
Macro pour envoyer une exception avec formattage.
#define ARCANE_FATAL(...)
Macro envoyant une exception FatalErrorException.
#define ARCANE_ALIGNAS_PACKED(value)
Macro pour garantir le compactage et l'alignement d'une classe sur value octets.
Classes, Types et macros pour gérer la concurrence.
#define ARCANE_REGISTER_APPLICATION_FACTORY(aclass, ainterface, aname)
Enregistre un service de fabrique pour la classe aclass.
Interval d'itération complexe.
Classe pour gérer le profiling d'une seule exécution d'une boucle.
Definition Profiling.h:93
Intervalle d'itération pour une boucle.
Informations d'exécution d'une boucle.
Interface d'un fonctor sur un interval d'itération multi-dimensionnel de dimension RankValue.
Interface d'un fonctor sur un interval d'itération.
virtual void executeFunctor(Integer begin, Integer size)=0
Exécute la méthode associée.
Interface d'un fonctor pour une tâche.
Implémentation d'une fabrique de tâches.
Interface d'une tâche concourante.
Fonctor sur un interval d'itération instancié via une lambda fonction.
void launchAndWait() override
Lance la tâche et bloque jusqu'à ce qu'elle se termine.
Lecteur des fichiers de maillage via la bibliothèque LIMA.
Definition Lima.cc:120
Caractéristiques d'un boucle 1D multi-thread.
Options d'exécution d'une boucle parallèle en multi-thread.
Integer grainSize() const
Taille d'un intervalle d'itération.
Int32 maxThread() const
Nombre maximal de threads autorisés.
static impl::ForLoopStatInfoList * _threadLocalForLoopInstance()
Definition Profiling.cc:194
static bool hasProfiling()
Indique si le profilage est actif.
Definition Profiling.h:173
Structure contenant les informations pour créer un service.
Implémentation déterministe de ParallelFor.
void operator()(tbb::blocked_range< Integer > &range) const
Opérateur pour un thread donné.
Exécuteur pour une boucle multi-dimension.
Exécuteur pour une boucle 1D.
std::vector< tbb::task_arena * > m_sub_arena_list
Tableau dont le i-ème élément contient la tbb::task_arena pour i thread.
Classe pour positionner TaskThreadInfo::taskIndex().
void executeParallelFor(const ComplexForLoopRanges< 1 > &loop_ranges, const ParallelLoopOptions &options, IMDRangeFunctor< 1 > *functor) final
Exécute une boucle 1D en concurrence.
Int32 currentTaskThreadIndex() const final
Implémentation de TaskFactory::currentTaskThreadIndex()
void initialize(Int32 nb_thread) override
void executeParallelFor(const ComplexForLoopRanges< 4 > &loop_ranges, const ParallelLoopOptions &options, IMDRangeFunctor< 4 > *functor) final
Exécute une boucle 4D en concurrence.
void executeParallelFor(const ComplexForLoopRanges< 3 > &loop_ranges, const ParallelLoopOptions &options, IMDRangeFunctor< 3 > *functor) final
Exécute une boucle 3D en concurrence.
void executeParallelFor(const ComplexForLoopRanges< 2 > &loop_ranges, const ParallelLoopOptions &options, IMDRangeFunctor< 2 > *functor) final
Exécute une boucle 2D en concurrence.
ITask * createRootTask(ITaskFunctor *f) override
Créé une tâche racine. L'implémentation doit recopier la valeur de f qui est soit un TaskFunctor,...
void printInfos(std::ostream &o) const final
Affiche les informations sur le runtime utilisé
void _executeMDParallelFor(const ComplexForLoopRanges< RankValue > &loop_ranges, IMDRangeFunctor< RankValue > *functor, const ParallelLoopOptions &options)
Exécution d'une boucle N-dimensions.
TaskThreadInfo * currentTaskThreadInfo() const
Instance de TaskThreadInfo associé au thread courant.
bool isActive() const final
Indique si l'implémentation est active.
Int32 currentTaskIndex() const final
Implémentation de TaskFactory::currentTaskIndex()
Int32 nbAllowedThread() const final
Nombre de threads utilisés au maximum pour gérer les tâches.
Contexte d'éxecution d'une tâche.
static IObservable * destroyThreadObservable()
Observable appelé lors de la destruction d'un thread pour une tâche.
static const ParallelLoopOptions & defaultParallelLoopOptions()
Valeurs par défaut d'exécution d'une boucle parallèle.
static IObservable * createThreadObservable()
Observable appelé lors de la création d'un thread pour une tâche.
static Integer verboseLevel()
Niveau de verbosité
static void setDefaultParallelLoopOptions(const ParallelLoopOptions &v)
Positionne les valeurs par défaut d'exécution d'une boucle parallèle.
Classe permettant de récupérer le temps passé entre l'appel au constructeur et au destructeur.
Definition Profiling.h:38
Allocateur mémoire avec alignement mémoire spécifique.
Exception lorsqu'une fonction n'est pas implémentée.
Integer toInteger(Real r)
Converti un Int64 en un Integer.
-*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
SimpleForLoopRanges< 1 > makeLoopRanges(Int32 n1)
Créé un intervalle d'itération [0,n1[.
Int32 Integer
Type représentant un entier.