Arcane  v3.15.3.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-2025 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-2025 */
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#include "arcane/utils/FixedArray.h"
25#include "arcane/utils/internal/TaskFactoryInternal.h"
26#include "arcane/utils/internal/DependencyInjection.h"
27
28#include "arcane/core/FactoryService.h"
29
30#include <new>
31#include <stack>
32
33// Il faut définir cette macro pour que la classe 'blocked_rangeNd' soit disponible
34
35#define TBB_PREVIEW_BLOCKED_RANGE_ND 1
36
37// la macro 'ARCANE_USE_ONETBB' est définie dans le CMakeLists.txt
38// si on compile avec la version OneTBB version 2021+
39// (https://github.com/oneapi-src/oneTBB.git)
40// A terme ce sera la seule version supportée par Arcane.
41
42#ifdef ARCANE_USE_ONETBB
43
44// Nécessaire pour avoir accès à task_scheduler_handle
45#define TBB_PREVIEW_WAITING_FOR_WORKERS 1
46#include <tbb/tbb.h>
47#include <oneapi/tbb/concurrent_set.h>
48#include <oneapi/tbb/global_control.h>
49
50#else // ARCANE_USE_ONETBB
51
52// NOTE GG: depuis mars 2019, la version 2018+ des TBB est obligatoire.
53#include <tbb/tbb.h>
54#include <tbb/blocked_rangeNd.h>
55/*
56 * Maintenant vérifie que la version est au moins 2018.
57 */
58#if (TBB_VERSION_MAJOR < 4) || (TBB_VERSION_MAJOR==4 && TBB_VERSION_MINOR<2)
59#define ARCANE_OLD_TBB
60#endif
61
62#if (TBB_VERSION_MAJOR<2018)
63#define ARCANE_OLD_TBB
64#endif
65
66#ifdef ARCANE_OLD_TBB
67#error "Your version of TBB is tool old. TBB 2018+ is required. Please disable TBB in configuration"
68#endif
69
70#include <thread>
71#include <mutex>
72
73#endif // ARCANE_USE_ONETBB
74
75/*---------------------------------------------------------------------------*/
76/*---------------------------------------------------------------------------*/
77
78namespace Arcane
79{
80
81class TBBTaskImplementation;
82
83// TODO: utiliser un pool mémoire spécifique pour gérer les
84// OneTBBTask pour optimiser les new/delete des instances de cette classe.
85// Auparavant avec les anciennes versions de TBB cela était géré avec
86// la méthode 'tbb::task::allocate_child()'.
87
88/*---------------------------------------------------------------------------*/
89/*---------------------------------------------------------------------------*/
90
91namespace
92{
93
94// Positif si on récupère les statistiques d'exécution
95bool isStatActive()
96{
98}
99
104class ScopedExecInfo
105{
106 public:
107
108 explicit ScopedExecInfo(const ForLoopRunInfo& run_info)
109 : m_run_info(run_info)
110 {
111 // Si run_info.execInfo() n'est pas nul, on l'utilise.
112 // Cela signifie que c'est l'appelant de qui va gérer les statistiques
113 // d'exécution. Sinon, on utilise \a m_stat_info si les statistiques
114 // d'exécution sont demandées.
115 ForLoopOneExecStat* ptr = run_info.execStat();
116 if (ptr){
117 m_stat_info_ptr = ptr;
118 m_use_own_run_info = false;
119 }
120 else
121 m_stat_info_ptr = isStatActive() ? &m_stat_info : nullptr;
122 }
123 ~ScopedExecInfo()
124 {
125#ifdef PRINT_STAT_INFO
126 if (m_stat_info_ptr){
127 bool is_valid = m_run_info.traceInfo().isValid();
128 if (!is_valid)
129 std::cout << "ADD_OWN_RUN_INFO nb_chunk=" << m_stat_info_ptr->nbChunk()
130 << " stack=" << platform::getStackTrace()
131 << "\n";
132 else
133 std::cout << "ADD_OWN_RUN_INFO nb_chunk=" << m_stat_info_ptr->nbChunk()
134 << " trace_name=" << m_run_info.traceInfo().traceInfo().name() << "\n";
135 }
136#endif
137 if (m_stat_info_ptr && m_use_own_run_info){
138 ProfilingRegistry::_threadLocalForLoopInstance()->merge(*m_stat_info_ptr,m_run_info.traceInfo());
139 }
140 }
141
142 public:
143
144 ForLoopOneExecStat* statInfo() const { return m_stat_info_ptr; }
145 bool isOwn() const { return m_use_own_run_info; }
146 private:
147
148 ForLoopOneExecStat m_stat_info;
149 ForLoopOneExecStat* m_stat_info_ptr = nullptr;
150 ForLoopRunInfo m_run_info;
152 bool m_use_own_run_info = true;
153};
154
155/*---------------------------------------------------------------------------*/
156/*---------------------------------------------------------------------------*/
157
158inline int _currentTaskTreadIndex()
159{
160 // NOTE: Avec OneTBB 2021, la valeur n'est plus '0' si on appelle cette méthode
161 // depuis un thread en dehors d'un task_arena. Avec la version 2021,
162 // la valeur est 65535.
163 // NOTE: Il semble que cela soit un bug de la 2021.3.
164 return tbb::this_task_arena::current_thread_index();
165}
166
167inline tbb::blocked_rangeNd<Int32,1>
168_toTBBRange(const ComplexForLoopRanges<1>& r)
169{
170 return {{r.lowerBound<0>(), r.upperBound<0>()}};
171}
172
173inline tbb::blocked_rangeNd<Int32,2>
174_toTBBRange(const ComplexForLoopRanges<2>& r)
175{
176 return {{r.lowerBound<0>(), r.upperBound<0>()},
177 {r.lowerBound<1>(), r.upperBound<1>()}};
178
179}
180
181inline tbb::blocked_rangeNd<Int32,3>
182_toTBBRange(const ComplexForLoopRanges<3>& r)
183{
184 return {{r.lowerBound<0>(), r.upperBound<0>()},
185 {r.lowerBound<1>(), r.upperBound<1>()},
186 {r.lowerBound<2>(), r.upperBound<2>()}};
187}
188
189inline tbb::blocked_rangeNd<Int32,4>
190_toTBBRange(const ComplexForLoopRanges<4>& r)
191{
192 return {{r.lowerBound<0>(), r.upperBound<0>()},
193 {r.lowerBound<1>(), r.upperBound<1>()},
194 {r.lowerBound<2>(), r.upperBound<2>()},
195 {r.lowerBound<3>(), r.upperBound<3>()}};
196}
197
198/*---------------------------------------------------------------------------*/
199/*---------------------------------------------------------------------------*/
200
201inline tbb::blocked_rangeNd<Int32,2>
202_toTBBRangeWithGrain(const tbb::blocked_rangeNd<Int32,2>& r,FixedArray<size_t,2> grain_sizes)
203{
204 return {{r.dim(0).begin(), r.dim(0).end(), grain_sizes[0]},
205 {r.dim(1).begin(), r.dim(1).end(), grain_sizes[1]}};
206}
207
208inline tbb::blocked_rangeNd<Int32,3>
209_toTBBRangeWithGrain(const tbb::blocked_rangeNd<Int32,3>& r,FixedArray<size_t,3> grain_sizes)
210{
211 return {{r.dim(0).begin(), r.dim(0).end(), grain_sizes[0]},
212 {r.dim(1).begin(), r.dim(1).end(), grain_sizes[1]},
213 {r.dim(2).begin(), r.dim(2).end(), grain_sizes[2]}};
214}
215
216inline tbb::blocked_rangeNd<Int32,4>
217_toTBBRangeWithGrain(const tbb::blocked_rangeNd<Int32,4>& r,FixedArray<size_t,4> grain_sizes)
218{
219 return {{r.dim(0).begin(), r.dim(0).end(), grain_sizes[0]},
220 {r.dim(1).begin(), r.dim(1).end(), grain_sizes[1]},
221 {r.dim(2).begin(), r.dim(2).end(), grain_sizes[2]},
222 {r.dim(3).begin(), r.dim(3).end(), grain_sizes[3]}};
223}
224
225/*---------------------------------------------------------------------------*/
226/*---------------------------------------------------------------------------*/
227
228inline ComplexForLoopRanges<2>
229_fromTBBRange(const tbb::blocked_rangeNd<Int32,2>& r)
230{
231 using BoundsType = ArrayBounds<MDDim2>;
232 using ArrayExtentType = typename BoundsType::ArrayExtentType;
233
234 BoundsType lower_bounds(ArrayExtentType(r.dim(0).begin(),r.dim(1).begin()));
235 auto s0 = static_cast<Int32>(r.dim(0).size());
236 auto s1 = static_cast<Int32>(r.dim(1).size());
237 BoundsType sizes(ArrayExtentType(s0,s1));
238 return { lower_bounds, sizes };
239}
240
241inline ComplexForLoopRanges<3>
242_fromTBBRange(const tbb::blocked_rangeNd<Int32,3> & r)
243{
244 using BoundsType = ArrayBounds<MDDim3>;
245 using ArrayExtentType = typename BoundsType::ArrayExtentType;
246
247 BoundsType lower_bounds(ArrayExtentType(r.dim(0).begin(),r.dim(1).begin(),r.dim(2).begin()));
248 auto s0 = static_cast<Int32>(r.dim(0).size());
249 auto s1 = static_cast<Int32>(r.dim(1).size());
250 auto s2 = static_cast<Int32>(r.dim(2).size());
251 BoundsType sizes(ArrayExtentType(s0,s1,s2));
252 return { lower_bounds, sizes };
253}
254
255inline ComplexForLoopRanges<4>
256_fromTBBRange(const tbb::blocked_rangeNd<Int32,4>& r)
257{
258 using BoundsType = ArrayBounds<MDDim4>;
259 using ArrayExtentType = typename BoundsType::ArrayExtentType;
260
261 BoundsType lower_bounds(ArrayExtentType(r.dim(0).begin(),r.dim(1).begin(),r.dim(2).begin(),r.dim(3).begin()));
262 auto s0 = static_cast<Int32>(r.dim(0).size());
263 auto s1 = static_cast<Int32>(r.dim(1).size());
264 auto s2 = static_cast<Int32>(r.dim(2).size());
265 auto s3 = static_cast<Int32>(r.dim(3).size());
266 BoundsType sizes(ArrayExtentType(s0,s1,s2,s3));
267 return { lower_bounds, sizes };
268}
269
270}
271
272/*---------------------------------------------------------------------------*/
273/*---------------------------------------------------------------------------*/
274
275#ifdef ARCANE_USE_ONETBB
276
277/*---------------------------------------------------------------------------*/
278/*---------------------------------------------------------------------------*/
279
280class OneTBBTaskFunctor
281{
282 public:
283 OneTBBTaskFunctor(ITaskFunctor* functor,ITask* task)
284 : m_functor(functor), m_task(task) {}
285 public:
286 void operator()() const
287 {
288 if (m_functor){
289 ITaskFunctor* tf = m_functor;
290 m_functor = 0;
291 TaskContext task_context(m_task);
292 //cerr << "FUNC=" << typeid(*tf).name();
293 tf->executeFunctor(task_context);
294 }
295 }
296 public:
297 mutable ITaskFunctor* m_functor;
298 ITask* m_task;
299};
300
301/*---------------------------------------------------------------------------*/
302/*---------------------------------------------------------------------------*/
303
304class OneTBBTask
305: public ITask
306{
307 public:
308 static const int FUNCTOR_CLASS_SIZE = 32;
309 public:
310 OneTBBTask(ITaskFunctor* f)
311 : m_functor(f)
312 {
313 m_functor = f->clone(functor_buf,FUNCTOR_CLASS_SIZE);
314 }
315 public:
316 OneTBBTaskFunctor taskFunctor() { return OneTBBTaskFunctor(m_functor,this); }
317 void launchAndWait() override;
318 void launchAndWait(ConstArrayView<ITask*> tasks) override;
319 protected:
320 virtual ITask* _createChildTask(ITaskFunctor* functor) override;
321 public:
322 ITaskFunctor* m_functor;
323 char functor_buf[FUNCTOR_CLASS_SIZE];
324};
325using TBBTask = OneTBBTask;
326
327/*---------------------------------------------------------------------------*/
328/*---------------------------------------------------------------------------*/
329
330#else // ARCANE_USE_ONETBB
331
332/*---------------------------------------------------------------------------*/
333/*---------------------------------------------------------------------------*/
334
336: public tbb::task
337, public ITask
338{
339 public:
340 static const int FUNCTOR_CLASS_SIZE = 32;
341 public:
343 : m_functor(f)
344 {
345 m_functor = f->clone(functor_buf,FUNCTOR_CLASS_SIZE);
346 }
347 public:
348 tbb::task* execute() override
349 {
350 if (m_functor){
351 ITaskFunctor* tf = m_functor;
352 m_functor = 0;
354 //cerr << "FUNC=" << typeid(*tf).name();
355 tf->executeFunctor(task_context);
356 }
357 return nullptr;
358 }
359 void launchAndWait() override;
361 protected:
362 ITask* _createChildTask(ITaskFunctor* functor) final;
363 public:
364 ITaskFunctor* m_functor;
365 char functor_buf[FUNCTOR_CLASS_SIZE];
366};
367using TBBTask = LegacyTBBTask;
368
369/*---------------------------------------------------------------------------*/
370/*---------------------------------------------------------------------------*/
371
372#endif // ARCANE_USE_ONETBB
373
374/*---------------------------------------------------------------------------*/
375/*---------------------------------------------------------------------------*/
376/*
377 * Ne pas utiliser l'observer locale au task_arena.
378 * Utiliser l'observer global au scheduler.
379 * Pour l'id, utiliser tbb::this_task_arena::current_thread_index().
380 */
382: public ITaskImplementation
383{
384 class Impl;
385 class ParallelForExecute;
386 template<int RankValue>
388
389 public:
390 // Pour des raisons de performance, s'aligne sur une ligne de cache
391 // et utilise un padding.
393 {
394 public:
395 TaskThreadInfo() : m_task_index(-1){}
396 public:
397 void setTaskIndex(Integer v) { m_task_index = v; }
398 Integer taskIndex() const { return m_task_index; }
399 private:
400 Integer m_task_index;
401 };
410 {
411 public:
413 : m_tti(tti), m_old_task_index(-1)
414 {
415 if (tti){
416 m_old_task_index = tti->taskIndex();
417 tti->setTaskIndex(task_index);
418 }
419 }
421 {
422 if (m_tti)
423 m_tti->setTaskIndex(m_old_task_index);
424 }
425 private:
426 TaskThreadInfo* m_tti;
427 Integer m_old_task_index;
428 };
429
430 public:
432 {
433 }
434 TBBTaskImplementation() = default;
435 ~TBBTaskImplementation() override;
436 public:
437 void build() {}
438 void initialize(Int32 nb_thread) override;
439 void terminate() override;
440
442 {
443#ifdef ARCANE_USE_ONETBB
444 OneTBBTask* t = new OneTBBTask(f);
445#else
446 LegacyTBBTask* t = new(tbb::task::allocate_root()) LegacyTBBTask(f);
447#endif
448 return t;
449 }
450
451 void executeParallelFor(Int32 begin,Int32 size,const ParallelLoopOptions& options,IRangeFunctor* f) final;
452 void executeParallelFor(Int32 begin,Int32 size,Integer grain_size,IRangeFunctor* f) final;
453 void executeParallelFor(Int32 begin,Int32 size,IRangeFunctor* f) final
454 {
455 executeParallelFor(begin,size,TaskFactory::defaultParallelLoopOptions(),f);
456 }
457 void executeParallelFor(const ParallelFor1DLoopInfo& loop_info) override;
458
483
484 bool isActive() const final
485 {
486 return m_is_active;
487 }
488
489 Int32 nbAllowedThread() const final;
490
492 {
493 return (nbAllowedThread() <= 1 ) ? 0 : _currentTaskTreadIndex();
494 }
495
496 Int32 currentTaskIndex() const final;
497
498 void printInfos(std::ostream& o) const final;
499
500 public:
501
508 TaskThreadInfo* currentTaskThreadInfo() const;
509
510 private:
511
512 bool m_is_active = false;
513 Impl* m_p = nullptr;
514
515 private:
516
517 template<int RankValue> void
519 IMDRangeFunctor<RankValue>* functor,
520 const ForLoopRunInfo& run_info);
521 void _executeParallelFor(const ParallelFor1DLoopInfo& loop_info);
522};
523
524/*---------------------------------------------------------------------------*/
525/*---------------------------------------------------------------------------*/
526
528{
530 : public tbb::task_scheduler_observer
531 {
532 public:
534 :
535#ifdef ARCANE_USE_ONETBB
536 tbb::task_scheduler_observer(p->m_main_arena),
537#endif
538 m_p(p)
539 {
540 }
541 void on_scheduler_entry(bool is_worker) override
542 {
543 m_p->notifyThreadCreated(is_worker);
544 }
545 void on_scheduler_exit(bool is_worker) override
546 {
547 m_p->notifyThreadDestroyed(is_worker);
548 }
550 };
551
552 public:
553 Impl() :
554 m_task_observer(this),
555 m_thread_task_infos(AlignedMemoryAllocator::CacheLine())
556 {
557#ifdef ARCANE_USE_ONETBB
558 m_nb_allowed_thread = tbb::info::default_concurrency();
559#else
560 m_nb_allowed_thread = tbb::task_scheduler_init::default_num_threads();
561#endif
562 _init();
563 }
564 Impl(Int32 nb_thread)
565 :
567 m_scheduler_init(nb_thread),
568#endif
569 m_main_arena(nb_thread),
570 m_task_observer(this),
571 m_thread_task_infos(AlignedMemoryAllocator::CacheLine())
572 {
573 m_nb_allowed_thread = nb_thread;
574 _init();
575 }
576 public:
577 Int32 nbAllowedThread() const { return m_nb_allowed_thread; }
578 TaskThreadInfo* threadTaskInfo(Integer index) { return &m_thread_task_infos[index]; }
579 private:
580 Int32 m_nb_allowed_thread;
581
582 public:
583 void terminate()
584 {
585 for( auto x : m_sub_arena_list ){
586 if (x)
587 x->terminate();
588 delete x;
589 }
590 m_sub_arena_list.clear();
591 m_main_arena.terminate();
592#ifdef ARCANE_USE_ONETBB
593 m_task_observer.observe(false);
594 oneapi::tbb::finalize(m_task_scheduler_handle);
595#else
596 m_scheduler_init.terminate();
597 m_task_observer.observe(false);
598#endif
599 }
600 public:
601 void notifyThreadCreated(bool is_worker)
602 {
603 std::thread::id my_thread_id = std::this_thread::get_id();
604
605#ifdef ARCANE_USE_ONETBB
606 // Avec OneTBB, cette méthode est appelée à chaque fois qu'on rentre
607 // dans notre 'task_arena'. Comme il ne faut appeler qu'une seule
608 // fois la méthode de notification on utilise un ensemble pour
609 // conserver la liste des threads déjà créés.
610 // NOTE: On ne peut pas utiliser cette méthode avec la version TBB historique
611 // (2018) car cette méthode 'contains' n'existe pas
612 if (m_constructed_thread_map.contains(my_thread_id))
613 return;
614 m_constructed_thread_map.insert(my_thread_id);
615#endif
616
617 {
619 std::ostringstream ostr;
620 ostr << "TBB: CREATE THREAD"
621 << " nb_allowed=" << m_nb_allowed_thread
622#ifdef ARCANE_USE_ONETBB
623 << " tbb_default_allowed=" << tbb::info::default_concurrency()
624#else
625 << " tbb_default_allowed=" << tbb::task_scheduler_init::default_num_threads()
626#endif
627 << " id=" << my_thread_id
628 << " arena_id=" << _currentTaskTreadIndex()
629 << " is_worker=" << is_worker
630 << "\n";
631 std::cout << ostr.str();
632 }
634 }
635 }
636
637 void notifyThreadDestroyed([[maybe_unused]] bool is_worker)
638 {
639#ifdef ARCANE_USE_ONETBB
640 // Avec OneTBB, cette méthode est appelée à chaque fois qu'on sort
641 // de l'arène principale. Du coup elle ne correspond pas vraiment à une
642 // destruction de thread. On ne fait donc rien pour cette notification
643 // TODO: Regarder comment on peut être notifié de la destruction effective
644 // du thread.
645#else
646 // Il faut toujours un verrou car on n'est pas certain que
647 // les méthodes appelées par l'observable soient thread-safe
648 // (et aussi TaskFactory::createThreadObservable() ne l'est pas)
649 std::scoped_lock sl(m_thread_created_mutex);
651 std::cout << "TBB: DESTROY THREAD"
652 << " id=" << std::this_thread::get_id()
653 << " arena_id=" << _currentTaskTreadIndex()
654 << " is_worker=" << is_worker
655 << '\n';
656 }
657 // TODO: jamais utilisé. Sera supprimé au passage à OneTBB.
658 TaskFactory::destroyThreadObservable()->notifyAllObservers();
659#endif
660 }
661 private:
662#ifdef ARCANE_USE_ONETBB
663#if TBB_VERSION_MAJOR>2021 || (TBB_VERSION_MAJOR==2021 && TBB_VERSION_MINOR>5)
664 oneapi::tbb::task_scheduler_handle m_task_scheduler_handle = oneapi::tbb::attach();
665#else
666 oneapi::tbb::task_scheduler_handle m_task_scheduler_handle = tbb::task_scheduler_handle::get();
667#endif
668#else
669 tbb::task_scheduler_init m_scheduler_init;
670#endif
671 public:
672 tbb::task_arena m_main_arena;
675 private:
676 TaskObserver m_task_observer;
677 std::mutex m_thread_created_mutex;
678 UniqueArray<TaskThreadInfo> m_thread_task_infos;
679#ifdef ARCANE_USE_ONETBB
680 tbb::concurrent_set<std::thread::id> m_constructed_thread_map;
681#endif
682 void _init()
683 {
685 std::cout << "TBB: TBBTaskImplementationInit nb_allowed_thread=" << m_nb_allowed_thread
686 << " id=" << std::this_thread::get_id()
687 << " version=" << TBB_VERSION_MAJOR << "." << TBB_VERSION_MINOR
688 << "\n";
689 }
690 m_thread_task_infos.resize(m_nb_allowed_thread);
691 m_task_observer.observe(true);
692 Integer max_arena_size = m_nb_allowed_thread;
693 // Limite artificiellement le nombre de tbb::task_arena
694 // pour éviter d'avoir trop d'objets alloués.
695 if (max_arena_size>512)
696 max_arena_size = 512;
697 if (max_arena_size<2)
698 max_arena_size = 2;
699 m_sub_arena_list.resize(max_arena_size);
700 m_sub_arena_list[0] = m_sub_arena_list[1] = nullptr;
701 for( Integer i=2; i<max_arena_size; ++i )
702 m_sub_arena_list[i] = new tbb::task_arena(i);
703 }
704};
705
706/*---------------------------------------------------------------------------*/
707/*---------------------------------------------------------------------------*/
712{
713 public:
715 : m_functor(f), m_stat_info(stat_info), m_nb_allowed_thread(nb_allowed_thread){}
716 public:
717
718 void operator()(tbb::blocked_range<Integer>& range) const
719 {
720#ifdef ARCANE_CHECK
721 if (TaskFactory::verboseLevel()>=3){
722 std::ostringstream o;
723 o << "TBB: INDEX=" << TaskFactory::currentTaskThreadIndex()
724 << " id=" << std::this_thread::get_id()
725 << " max_allowed=" << m_nb_allowed_thread
726 << " range_begin=" << range.begin() << " range_size=" << range.size()
727 << "\n";
728 std::cout << o.str();
729 std::cout.flush();
730 }
731
733 if (tbb_index<0 || tbb_index>=m_nb_allowed_thread)
734 ARCANE_FATAL("Invalid index for thread idx={0} valid_interval=[0..{1}[",
735 tbb_index,m_nb_allowed_thread);
736#endif
737
738 if (m_stat_info)
739 m_stat_info->incrementNbChunk();
740 m_functor->executeFunctor(range.begin(),CheckedConvert::toInteger(range.size()));
741 }
742
743 private:
744 IRangeFunctor* m_functor;
745 ForLoopOneExecStat* m_stat_info = nullptr;
746 Int32 m_nb_allowed_thread;
747};
748
749/*---------------------------------------------------------------------------*/
750/*---------------------------------------------------------------------------*/
754template<int RankValue>
756{
757 public:
758
760 : m_functor(f), m_stat_info(stat_info), m_nb_allowed_thread(nb_allowed_thread){}
761
762 public:
763
764 void operator()(tbb::blocked_rangeNd<Int32,RankValue>& range) const
765 {
766#ifdef ARCANE_CHECK
767 if (TaskFactory::verboseLevel()>=3){
768 std::ostringstream o;
769 o << "TBB: INDEX=" << TaskFactory::currentTaskThreadIndex()
770 << " id=" << std::this_thread::get_id()
771 << " max_allowed=" << m_nb_allowed_thread
772 << " MDFor ";
773 for( Int32 i=0; i<RankValue; ++i ){
774 Int32 r0 = static_cast<Int32>(range.dim(i).begin());
775 Int32 r1 = static_cast<Int32>(range.dim(i).size());
776 o << " range" << i << " (begin=" << r0 << " size=" << r1 << ")";
777 }
778 o << "\n";
779 std::cout << o.str();
780 std::cout.flush();
781 }
782
784 if (tbb_index<0 || tbb_index>=m_nb_allowed_thread)
785 ARCANE_FATAL("Invalid index for thread idx={0} valid_interval=[0..{1}[",
786 tbb_index,m_nb_allowed_thread);
787#endif
788
789 if (m_stat_info)
790 m_stat_info->incrementNbChunk();
791 m_functor->executeFunctor(_fromTBBRange(range));
792 }
793
794 private:
795
796 IMDRangeFunctor<RankValue>* m_functor = nullptr;
797 ForLoopOneExecStat* m_stat_info = nullptr;
798 Int32 m_nb_allowed_thread;
799};
800
801/*---------------------------------------------------------------------------*/
802/*---------------------------------------------------------------------------*/
821{
822 public:
824 Integer begin_index,Integer size,Integer grain_size,Integer nb_thread)
825 : m_impl(impl), m_tbb_for(tbb_for), m_nb_thread(nb_thread), m_begin_index(begin_index), m_size(size),
826 m_grain_size(grain_size), m_nb_block(0), m_block_size(0), m_nb_block_per_thread(0)
827 {
828 if (m_nb_thread<1)
829 m_nb_thread = 1;
830
831 if (m_grain_size>0){
832 m_block_size = m_grain_size;
833 if (m_block_size>0){
834 m_nb_block = m_size / m_block_size;
835 if ((m_size % m_block_size)!=0)
836 ++m_nb_block;
837 }
838 else
839 m_nb_block = 1;
840 m_nb_block_per_thread = m_nb_block / m_nb_thread;
841 if ((m_nb_block % m_nb_thread) != 0)
842 ++m_nb_block_per_thread;
843 }
844 else{
845 if (m_nb_block<1)
846 m_nb_block = m_nb_thread;
847 m_block_size = m_size / m_nb_block;
848 m_nb_block_per_thread = 1;
849 }
850 if (TaskFactory::verboseLevel()>=2){
851 std::cout << "TBBDeterministicParallelFor: BEGIN=" << m_begin_index << " size=" << m_size
852 << " grain_size=" << m_grain_size
853 << " nb_block=" << m_nb_block << " nb_thread=" << m_nb_thread
854 << " nb_block_per_thread=" << m_nb_block_per_thread
855 << " block_size=" << m_block_size
856 << " block_size*nb_block=" << m_block_size*m_nb_block << '\n';
857 }
858 }
859 public:
860
867 void operator()(tbb::blocked_range<Integer>& range) const
868 {
869 Integer nb_iter = static_cast<Integer>(range.size());
870 for( Integer i=0; i<nb_iter; ++i ){
871 Integer task_id = range.begin() + i;
872 for ( Integer k=0, kn=m_nb_block_per_thread; k<kn; ++k ){
873 Integer block_id = task_id + (k * m_nb_thread);
874 if (block_id<m_nb_block)
875 _doBlock(task_id,block_id);
876 }
877 }
878 }
879
880 void _doBlock(Integer task_id,Integer block_id) const
881 {
882 TBBTaskImplementation::TaskInfoLockGuard guard(m_impl->currentTaskThreadInfo(),task_id);
883
884 Integer iter_begin = block_id * m_block_size;
885 Integer iter_size = m_block_size;
886 if ((block_id+1)==m_nb_block){
887 // Pour le dernier bloc, la taille est le nombre d'éléments restants
888 iter_size = m_size - iter_begin;
889 }
890 iter_begin += m_begin_index;
891#ifdef ARCANE_CHECK
892 if (TaskFactory::verboseLevel()>=3){
893 std::ostringstream o;
894 o << "TBB: DoBlock: BLOCK task_id=" << task_id << " block_id=" << block_id
895 << " iter_begin=" << iter_begin << " iter_size=" << iter_size << '\n';
896 std::cout << o.str();
897 std::cout.flush();
898 }
899#endif
900 if (iter_size>0){
901 auto r = tbb::blocked_range<int>(iter_begin,iter_begin + iter_size);
902 m_tbb_for(r);
903 }
904 }
905
906 private:
907
908 TBBTaskImplementation* m_impl;
909 const TBBParallelFor& m_tbb_for;
910 Integer m_nb_thread;
911 Integer m_begin_index;
912 Integer m_size;
913 Integer m_grain_size;
914 Integer m_nb_block;
915 Integer m_block_size;
916 Integer m_nb_block_per_thread;
917};
918
919/*---------------------------------------------------------------------------*/
920/*---------------------------------------------------------------------------*/
921
923{
924 public:
925
927 Integer begin,Integer size,IRangeFunctor* f,ForLoopOneExecStat* stat_info)
928 : m_impl(impl), m_begin(begin), m_size(size), m_functor(f), m_options(options), m_stat_info(stat_info){}
929
930 public:
931
932 void operator()() const
933 {
934 Integer nb_thread = m_options.maxThread();
935 TBBParallelFor pf(m_functor,nb_thread,m_stat_info);
936 Integer gsize = m_options.grainSize();
937 tbb::blocked_range<Integer> range(m_begin,m_begin+m_size);
938 if (TaskFactory::verboseLevel()>=1)
939 std::cout << "TBB: TBBTaskImplementationInit ParallelForExecute begin=" << m_begin
940 << " size=" << m_size << " gsize=" << gsize
941 << " partitioner=" << (int)m_options.partitioner()
942 << " nb_thread=" << nb_thread
943 << " has_stat_info=" << (m_stat_info!=nullptr)
944 << '\n';
945
946 if (gsize>0)
947 range = tbb::blocked_range<Integer>(m_begin,m_begin+m_size,gsize);
948
949 if (m_options.partitioner()==ParallelLoopOptions::Partitioner::Static){
950 tbb::parallel_for(range,pf,tbb::static_partitioner());
951 }
952 else if (m_options.partitioner()==ParallelLoopOptions::Partitioner::Deterministic){
953 tbb::blocked_range<Integer> range2(0,nb_thread,1);
954 TBBDeterministicParallelFor dpf(m_impl,pf,m_begin,m_size,gsize,nb_thread);
955 tbb::parallel_for(range2,dpf);
956 }
957 else
958 tbb::parallel_for(range,pf);
959 }
960 private:
961 TBBTaskImplementation* m_impl = nullptr;
962 Integer m_begin;
963 Integer m_size;
964 IRangeFunctor* m_functor = nullptr;
965 ParallelLoopOptions m_options;
966 ForLoopOneExecStat* m_stat_info = nullptr;
967};
968
969/*---------------------------------------------------------------------------*/
970/*---------------------------------------------------------------------------*/
971
972template<int RankValue>
974{
975 public:
976
978 const ParallelLoopOptions& options,
981 : m_impl(impl)
982 , m_tbb_range(_toTBBRange(range))
983 , m_functor(f)
984 , m_options(options)
985 , m_stat_info(stat_info)
986 {
987 // On ne peut pas modifier les valeurs d'une instance de tbb::blocked_rangeNd.
988 // Il faut donc en reconstruire une complètement.
990 Int32 gsize = m_options.grainSize();
991 if (gsize>0){
992 // Si la taille du grain est différent zéro, il faut la répartir
993 // sur l'ensemble des dimensions. On commence par la dernière.
994 // TODO: regarder pourquoi dans certains cas les performances sont
995 // inférieures à celles qu'on obtient en utilisant un partitionneur
996 // statique.
997 constexpr bool is_verbose = false;
998 std::array<Int32,RankValue> range_extents = range.extents().asStdArray();
999 double ratio = static_cast<double>(gsize) / static_cast<double>(range.nbElement());
1000 if constexpr (is_verbose){
1001 std::cout << "GSIZE=" << gsize << " rank=" << RankValue << " ratio=" << ratio;
1002 for(Int32 i=0; i<RankValue; ++i )
1003 std::cout << " range" << i << "=" << range_extents[i];
1004 std::cout << "\n";
1005 }
1006 Int32 index = RankValue - 1;
1007 Int32 remaining_grain = gsize;
1008 for( ; index>=0; --index ){
1009 Int32 current = range_extents[index];
1010 if constexpr (is_verbose)
1011 std::cout << "Check index=" << index << " remaining=" << remaining_grain << " current=" << current << "\n";
1012 if (remaining_grain>current){
1013 all_grain_sizes[index] = current;
1014 remaining_grain /= current;
1015 }
1016 else{
1018 break;
1019 }
1020 }
1021 for( Int32 i=0; i<index; ++i )
1022 all_grain_sizes[i] = 1;
1023 if constexpr (is_verbose){
1024 for(Int32 i=0; i<RankValue; ++i )
1025 std::cout << " grain" << i << "=" << all_grain_sizes[i];
1026 std::cout << "\n";
1027 }
1028 m_tbb_range = _toTBBRangeWithGrain(m_tbb_range,all_grain_sizes);
1029 }
1030 }
1031
1032 public:
1033
1034 void operator()() const
1035 {
1036 Integer nb_thread = m_options.maxThread();
1037 TBBMDParallelFor<RankValue> pf(m_functor,nb_thread,m_stat_info);
1038
1039 if (m_options.partitioner()==ParallelLoopOptions::Partitioner::Static){
1040 tbb::parallel_for(m_tbb_range,pf,tbb::static_partitioner());
1041 }
1042 else if (m_options.partitioner()==ParallelLoopOptions::Partitioner::Deterministic){
1043 // TODO: implémenter le mode déterministe
1044 ARCANE_THROW(NotImplementedException,"ParallelLoopOptions::Partitioner::Deterministic for multi-dimensionnal loops");
1045 //tbb::blocked_range<Integer> range2(0,nb_thread,1);
1046 //TBBDeterministicParallelFor dpf(m_impl,pf,m_begin,m_size,gsize,nb_thread);
1047 //tbb::parallel_for(range2,dpf);
1048 }
1049 else{
1050 tbb::parallel_for(m_tbb_range,pf);
1051 }
1052 }
1053 private:
1054 TBBTaskImplementation* m_impl = nullptr;
1055 tbb::blocked_rangeNd<Int32,RankValue> m_tbb_range;
1056 IMDRangeFunctor<RankValue>* m_functor = nullptr;
1057 ParallelLoopOptions m_options;
1058 ForLoopOneExecStat* m_stat_info = nullptr;
1059};
1060
1061/*---------------------------------------------------------------------------*/
1062/*---------------------------------------------------------------------------*/
1063
1064TBBTaskImplementation::
1065~TBBTaskImplementation()
1066{
1067 delete m_p;
1068}
1069
1070/*---------------------------------------------------------------------------*/
1071/*---------------------------------------------------------------------------*/
1072
1074initialize(Int32 nb_thread)
1075{
1076 if (nb_thread<0)
1077 nb_thread = 0;
1078 m_is_active = (nb_thread!=1);
1079 if (nb_thread!=0)
1080 m_p = new Impl(nb_thread);
1081 else
1082 m_p = new Impl();
1084 opts.setMaxThread(nbAllowedThread());
1086}
1087
1088/*---------------------------------------------------------------------------*/
1089/*---------------------------------------------------------------------------*/
1090
1092terminate()
1093{
1094 m_p->terminate();
1095}
1096
1097/*---------------------------------------------------------------------------*/
1098/*---------------------------------------------------------------------------*/
1099
1101nbAllowedThread() const
1102{
1103 return m_p->nbAllowedThread();
1104}
1105
1106/*---------------------------------------------------------------------------*/
1107/*---------------------------------------------------------------------------*/
1108
1110printInfos(std::ostream& o) const
1111{
1112#ifdef ARCANE_USE_ONETBB
1113 o << "OneTBBTaskImplementation"
1114 << " version=" << TBB_VERSION_STRING
1115 << " interface=" << TBB_INTERFACE_VERSION
1116 << " runtime_interface=" << TBB_runtime_interface_version();
1117#else
1118 o << "TBBTaskImplementation"
1119 << " version=" << TBB_VERSION_MAJOR << "." << TBB_VERSION_MINOR
1120 << " interface=" << TBB_INTERFACE_VERSION;
1121#endif
1122}
1123
1124/*---------------------------------------------------------------------------*/
1125/*---------------------------------------------------------------------------*/
1126
1127void TBBTaskImplementation::
1128_executeParallelFor(const ParallelFor1DLoopInfo& loop_info)
1129{
1130 ScopedExecInfo sei(loop_info.runInfo());
1131 ForLoopOneExecStat* stat_info = sei.statInfo();
1133
1134 Int32 begin = loop_info.beginIndex();
1135 Int32 size = loop_info.size();
1136 ParallelLoopOptions options = loop_info.runInfo().options().value_or(TaskFactory::defaultParallelLoopOptions());
1137 IRangeFunctor* f = loop_info.functor();
1138
1139 Integer max_thread = options.maxThread();
1140 Integer nb_allowed_thread = m_p->nbAllowedThread();
1141 if (max_thread<0)
1143
1145 std::cout << "TBB: TBBTaskImplementation executeParallelFor begin=" << begin
1146 << " size=" << size << " max_thread=" << max_thread
1147 << " grain_size=" << options.grainSize()
1148 << " nb_allowed=" << nb_allowed_thread << '\n';
1149
1150 // En exécution séquentielle, appelle directement la méthode \a f.
1151 if (max_thread==1 || max_thread==0){
1152 f->executeFunctor(begin,size);
1153 return;
1154 }
1155
1156 // Remplace les valeurs non initialisées de \a options par celles de \a m_default_loop_options
1157 ParallelLoopOptions true_options(options);
1158 true_options.mergeUnsetValues(TaskFactory::defaultParallelLoopOptions());
1159 true_options.setMaxThread(max_thread);
1160
1161 ParallelForExecute pfe(this,true_options,begin,size,f,stat_info);
1162
1163 tbb::task_arena* used_arena = nullptr;
1164 if (max_thread<nb_allowed_thread && max_thread<m_p->m_sub_arena_list.size())
1165 used_arena = m_p->m_sub_arena_list[max_thread];
1166 if (!used_arena)
1167 used_arena = &(m_p->m_main_arena);
1168 used_arena->execute(pfe);
1169}
1170
1171/*---------------------------------------------------------------------------*/
1172/*---------------------------------------------------------------------------*/
1173
1174void TBBTaskImplementation::
1175executeParallelFor(const ParallelFor1DLoopInfo& loop_info)
1176{
1177 _executeParallelFor(loop_info);
1178}
1179
1180/*---------------------------------------------------------------------------*/
1181/*---------------------------------------------------------------------------*/
1188template<int RankValue> void TBBTaskImplementation::
1191 const ForLoopRunInfo& run_info)
1192{
1193 ParallelLoopOptions options;
1194 if (run_info.options().has_value())
1195 options = run_info.options().value();
1196
1197 ScopedExecInfo sei(run_info);
1198 ForLoopOneExecStat* stat_info = sei.statInfo();
1199 impl::ScopedStatLoop scoped_loop(sei.isOwn() ? stat_info : nullptr);
1200
1201 if (TaskFactory::verboseLevel()>=1){
1202 std::cout << "TBB: TBBTaskImplementation executeMDParallelFor nb_dim=" << RankValue
1203 << " nb_element=" << loop_ranges.nbElement()
1204 << " grain_size=" << options.grainSize()
1205 << " name=" << run_info.traceInfo().traceInfo()
1206 << " has_stat_info=" << (stat_info!=nullptr)
1207 << '\n';
1208 }
1209
1210 Integer max_thread = options.maxThread();
1211 // En exécution séquentielle, appelle directement la méthode \a f.
1212 if (max_thread==1 || max_thread==0){
1213 functor->executeFunctor(loop_ranges);
1214 return;
1215 }
1216
1217 // Remplace les valeurs non initialisées de \a options par celles de \a m_default_loop_options
1220
1221 Integer nb_allowed_thread = m_p->nbAllowedThread();
1222 if (max_thread<0)
1224 tbb::task_arena* used_arena = nullptr;
1227 if (!used_arena)
1228 used_arena = &(m_p->m_main_arena);
1229
1230 // Pour l'instant pour la dimension 1, utilise le 'ParallelForExecute' historique
1231 if constexpr (RankValue==1){
1233 auto x1 = [&](Integer begin,Integer size)
1234 {
1235 functor->executeFunctor(makeLoopRanges(ForLoopRange(begin,size)));
1236 //functor->executeFunctor(ComplexForLoopRanges<1>(begin,size));
1237 };
1238 LambdaRangeFunctorT<decltype(x1)> functor_1d(x1);
1239 Integer begin1 = CheckedConvert::toInteger(range_1d.dim(0).begin());
1240 Integer size1 = CheckedConvert::toInteger(range_1d.dim(0).size());
1242 used_arena->execute(pfe);
1243 }
1244 else{
1246 used_arena->execute(pfe);
1247 }
1248}
1249
1250/*---------------------------------------------------------------------------*/
1251/*---------------------------------------------------------------------------*/
1252
1253void TBBTaskImplementation::
1254executeParallelFor(Integer begin,Integer size,Integer grain_size,IRangeFunctor* f)
1255{
1257 opts.setGrainSize(grain_size);
1259 executeParallelFor(ParallelFor1DLoopInfo(begin,size,f,run_info));
1260}
1261
1262/*---------------------------------------------------------------------------*/
1263/*---------------------------------------------------------------------------*/
1264
1265void TBBTaskImplementation::
1266executeParallelFor(Integer begin,Integer size,const ParallelLoopOptions& options,IRangeFunctor* f)
1267{
1268 executeParallelFor(ParallelFor1DLoopInfo(begin,size,f,ForLoopRunInfo(options)));
1269}
1270
1271/*---------------------------------------------------------------------------*/
1272/*---------------------------------------------------------------------------*/
1273
1276{
1278 if (thread_id>=0)
1279 return m_p->threadTaskInfo(thread_id);
1280 return nullptr;
1281}
1282
1283/*---------------------------------------------------------------------------*/
1284/*---------------------------------------------------------------------------*/
1285
1287currentTaskIndex() const
1288{
1290#ifdef ARCANE_USE_ONETBB
1291 if (thread_id<0 || thread_id>=m_p->nbAllowedThread())
1292 return 0;
1293#endif
1295 if (tti){
1296 Int32 task_index = tti->taskIndex();
1297 if (task_index>=0)
1298 return task_index;
1299 }
1300 return thread_id;
1301}
1302
1303/*---------------------------------------------------------------------------*/
1304/*---------------------------------------------------------------------------*/
1305
1306#ifdef ARCANE_USE_ONETBB
1307
1308/*---------------------------------------------------------------------------*/
1309/*---------------------------------------------------------------------------*/
1310
1311void OneTBBTask::
1312launchAndWait()
1313{
1314 tbb::task_group task_group;
1315 task_group.run(taskFunctor());
1316 task_group.wait();
1317 delete this;
1318}
1319
1320/*---------------------------------------------------------------------------*/
1321/*---------------------------------------------------------------------------*/
1322
1323void OneTBBTask::
1324launchAndWait(ConstArrayView<ITask*> tasks)
1325{
1326 tbb::task_group task_group;
1327 Integer n = tasks.size();
1328 if (n==0)
1329 return;
1330
1331 //set_ref_count(n+1);
1332 for( Integer i=0; i<n; ++i ){
1333 OneTBBTask* t = static_cast<OneTBBTask*>(tasks[i]);
1334 task_group.run(t->taskFunctor());
1335 }
1336 task_group.wait();
1337 for( Integer i=0; i<n; ++i ){
1338 OneTBBTask* t = static_cast<OneTBBTask*>(tasks[i]);
1339 delete t;
1340 }
1341}
1342
1343/*---------------------------------------------------------------------------*/
1344/*---------------------------------------------------------------------------*/
1345
1346ITask* OneTBBTask::
1347_createChildTask(ITaskFunctor* functor)
1348{
1349 OneTBBTask* t = new OneTBBTask(functor);
1350 return t;
1351}
1352
1353/*---------------------------------------------------------------------------*/
1354/*---------------------------------------------------------------------------*/
1355
1356#else // ARCANE_USE_ONETBB
1357
1358/*---------------------------------------------------------------------------*/
1359/*---------------------------------------------------------------------------*/
1360
1363{
1364 task::spawn_root_and_wait(*this);
1365}
1366
1367/*---------------------------------------------------------------------------*/
1368/*---------------------------------------------------------------------------*/
1369
1372{
1373 Integer n = tasks.size();
1374 if (n==0)
1375 return;
1376
1377 set_ref_count(n+1);
1378 for( Integer i=0; i<n-1; ++i ){
1379 TBBTask* t = static_cast<TBBTask*>(tasks[i]);
1380 spawn(*t);
1381 }
1382 spawn_and_wait_for_all(*static_cast<TBBTask*>(tasks[n-1]));
1383}
1384
1385/*---------------------------------------------------------------------------*/
1386/*---------------------------------------------------------------------------*/
1387
1388ITask* LegacyTBBTask::
1389_createChildTask(ITaskFunctor* functor)
1390{
1391 TBBTask* t = new(allocate_child()) TBBTask(functor);
1392 return t;
1393}
1394
1395/*---------------------------------------------------------------------------*/
1396/*---------------------------------------------------------------------------*/
1397
1398#endif // ARCANE_USE_ONETBB
1399
1400/*---------------------------------------------------------------------------*/
1401/*---------------------------------------------------------------------------*/
1402
1403// TODO: a supprimer maintenant qu'on utilise 'DependencyInjection'
1404ARCANE_REGISTER_APPLICATION_FACTORY(TBBTaskImplementation,ITaskImplementation,
1405 TBBTaskImplementation);
1406
1407ARCANE_DI_REGISTER_PROVIDER(TBBTaskImplementation,
1408 DependencyInjection::ProviderProperty("TBBTaskImplementation"),
1409 ARCANE_DI_INTERFACES(ITaskImplementation),
1410 ARCANE_DI_EMPTY_CONSTRUCTOR());
1411
1412/*---------------------------------------------------------------------------*/
1413/*---------------------------------------------------------------------------*/
1414
1415} // End namespace Arcane
1416
1417/*---------------------------------------------------------------------------*/
1418/*---------------------------------------------------------------------------*/
#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:149
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.
UniqueArray< 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().
Int32 currentTaskThreadIndex() const final
Implémentation de TaskFactory::currentTaskThreadIndex()
void initialize(Int32 nb_thread) override
void executeParallelFor(const ComplexForLoopRanges< 1 > &loop_ranges, const ForLoopRunInfo &run_info, IMDRangeFunctor< 1 > *functor) final
Exécute une boucle 1D 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 executeParallelFor(const ComplexForLoopRanges< 3 > &loop_ranges, const ForLoopRunInfo &run_info, IMDRangeFunctor< 3 > *functor) final
Exécute une boucle 3D en concurrence.
void executeParallelFor(const ComplexForLoopRanges< 4 > &loop_ranges, const ForLoopRunInfo &run_info, IMDRangeFunctor< 4 > *functor) final
Exécute une boucle 4D en concurrence.
TaskThreadInfo * currentTaskThreadInfo() const
Instance de TaskThreadInfo associé au thread courant.
bool isActive() const final
Indique si l'implémentation est active.
void _executeMDParallelFor(const ComplexForLoopRanges< RankValue > &loop_ranges, IMDRangeFunctor< RankValue > *functor, const ForLoopRunInfo &run_info)
Exécution d'une boucle N-dimensions.
Int32 currentTaskIndex() const final
Implémentation de TaskFactory::currentTaskIndex()
void executeParallelFor(const ComplexForLoopRanges< 2 > &loop_ranges, const ForLoopRunInfo &run_info, IMDRangeFunctor< 2 > *functor) final
Exécute une boucle 2D en concurrence.
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 void notifyThreadCreated()
Notifie tous les observateurs de création de thread.
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 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.