18 #ifndef LIBMESH_THREADS_PTHREAD_H 19 #define LIBMESH_THREADS_PTHREAD_H 23 #ifndef LIBMESH_SQUASH_HEADER_WARNING 24 # warning "This file is designed to be included through libmesh/threads.h" 27 #ifdef LIBMESH_HAVE_PTHREAD 30 #ifdef LIBMESH_HAVE_CXX11_THREAD 43 # include <libkern/OSAtomic.h> 48 #ifdef LIBMESH_HAVE_CXX11_THREAD 49 # define LIBMESH_TLS_TYPE(type) thread_local type 50 # define LIBMESH_TLS_REF(value) (value) 51 #else // Maybe support gcc __thread eventually? 52 # define LIBMESH_TLS_TYPE(type) type 53 # define LIBMESH_TLS_REF(value) (value) 63 #ifdef LIBMESH_HAVE_CXX11_THREAD 67 typedef std::thread
Thread;
74 typedef NonConcurrentThread
Thread;
76 #endif // LIBMESH_HAVE_CXX11_THREAD 184 pthread_mutexattr_init(&
attr);
185 pthread_mutexattr_settype(&
attr, PTHREAD_MUTEX_RECURSIVE);
213 template <
typename Range>
217 return min > 0 ? cast_int<unsigned int>(min) : 1;
220 template <
typename Range,
typename Body>
228 template <
typename Range,
typename Body>
233 Body & body = *range_body->
body;
234 Range & range = *range_body->
range;
268 template <
typename Range,
typename Body>
270 void parallel_for (
const Range & range,
const Body & body)
281 DisablePerfLogInScope disable_perf;
285 std::vector<std::unique_ptr<Range>> ranges(
n_threads);
286 std::vector<RangeBody<const Range, const Body>> range_bodies(
n_threads);
287 std::vector<pthread_t> threads(
n_threads);
290 std::size_t range_size = range.size() /
n_threads;
292 typename Range::const_iterator current_beginning = range.begin();
296 std::size_t this_range_size = range_size;
299 this_range_size += range.size() %
n_threads;
301 ranges[i] = std::make_unique<Range>(range, current_beginning, current_beginning + this_range_size);
303 current_beginning = current_beginning + this_range_size;
309 range_bodies[i].range = ranges[i].get();
310 range_bodies[i].body = &body;
317 #ifdef LIBMESH_HAVE_OPENMP 318 #pragma omp parallel for schedule (static) 320 for (
int i=0; i<static_cast<int>(
n_threads); i++)
322 #if !LIBMESH_HAVE_OPENMP 323 pthread_create(&threads[i],
nullptr, &run_body<Range, Body>, (
void *)&range_bodies[i]);
325 run_body<Range, Body>((
void *)&range_bodies[i]);
329 #if !LIBMESH_HAVE_OPENMP 338 for (
int i=0; i<static_cast<int>(
n_threads); i++)
339 pthread_join(threads[i],
nullptr);
347 template <
typename Range,
typename Body,
typename Partitioner>
349 void parallel_for (
const Range & range,
const Body & body,
const Partitioner &)
358 template <
typename Range,
typename Body>
371 DisablePerfLogInScope disable_perf;
375 std::vector<std::unique_ptr<Range>> ranges(
n_threads);
376 std::vector<std::unique_ptr<Body>> managed_bodies(
n_threads);
378 std::vector<RangeBody<Range, Body>> range_bodies(
n_threads);
383 managed_bodies[i] = std::make_unique<Body>(body,
Threads::split());
389 bodies[i] = managed_bodies[i].
get();
392 std::size_t range_size = range.size() /
n_threads;
394 typename Range::const_iterator current_beginning = range.begin();
398 std::size_t this_range_size = range_size;
401 this_range_size += range.size() %
n_threads;
403 ranges[i] = std::make_unique<Range>(range, current_beginning, current_beginning + this_range_size);
405 current_beginning = current_beginning + this_range_size;
411 range_bodies[i].range = ranges[i].get();
412 range_bodies[i].body = bodies[i];
416 std::vector<pthread_t> threads(
n_threads);
421 #ifdef LIBMESH_HAVE_OPENMP 422 #pragma omp parallel for schedule (static) 430 for (
int i=0; i<static_cast<int>(
n_threads); i++)
432 #if !LIBMESH_HAVE_OPENMP 433 pthread_create(&threads[i],
nullptr, &run_body<Range, Body>, (
void *)&range_bodies[i]);
435 run_body<Range, Body>((
void *)&range_bodies[i]);
439 #if !LIBMESH_HAVE_OPENMP 442 pthread_join(threads[i],
nullptr);
446 for (
unsigned int i=
n_threads-1; i != 0; i--)
447 bodies[i-1]->join(*bodies[i]);
454 template <
typename Range,
typename Body,
typename Partitioner>
456 void parallel_reduce (
const Range & range, Body & body,
const Partitioner &)
466 template <
typename T>
471 operator T () {
return val; }
475 spin_mutex::scoped_lock lock(
smutex);
482 spin_mutex::scoped_lock lock(
smutex);
490 spin_mutex::scoped_lock lock(
smutex);
497 spin_mutex::scoped_lock lock(
smutex);
504 spin_mutex::scoped_lock lock(
smutex);
511 spin_mutex::scoped_lock lock(
smutex);
518 spin_mutex::scoped_lock lock(
smutex);
525 spin_mutex::scoped_lock lock(
smutex);
539 #endif // #ifdef LIBMESH_HAVE_PTHREAD 541 #endif // LIBMESH_SQUASH_HEADER_WARNING 543 #endif // LIBMESH_THREADS_PTHREAD_H
scoped_lock(recursive_mutex &in_rmutex)
NonConcurrentThread Thread
Use the non-concurrent placeholder.
void acquire(recursive_mutex &in_rmutex)
void initialize(int=automatic)
void parallel_for(const Range &range, const Body &body)
Execute the provided function object in parallel on the specified range.
void * run_body(void *args)
The libMesh namespace provides an interface to certain functionality in the library.
bool in_threads
A boolean which is true iff we are in a Threads:: function It may be useful to assert(!Threadsin_thre...
tbb::task_scheduler_init task_scheduler_init
Scheduler to manage threads.
scoped_lock(spin_mutex &in_smutex)
void acquire(spin_mutex &in_smutex)
tbb::spin_mutex spin_mutex
Spin mutex.
tbb::split split
Dummy "splitting object" used to distinguish splitting constructors from copy constructors.
Defines atomic operations which can only be executed on a single thread at a time.
atomic< T > & operator=(const atomic< T > &value)
void parallel_reduce(const Range &range, Body &body)
Execute the provided reduction operation in parallel on the specified range.
static const int automatic
unsigned int num_pthreads(Range &range)
task_scheduler_init(int=automatic)