threads.h
Go to the documentation of this file.00001 // The libMesh Finite Element Library. 00002 // Copyright (C) 2002-2012 Benjamin S. Kirk, John W. Peterson, Roy H. Stogner 00003 00004 // This library is free software; you can redistribute it and/or 00005 // modify it under the terms of the GNU Lesser General Public 00006 // License as published by the Free Software Foundation; either 00007 // version 2.1 of the License, or (at your option) any later version. 00008 00009 // This library is distributed in the hope that it will be useful, 00010 // but WITHOUT ANY WARRANTY; without even the implied warranty of 00011 // MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU 00012 // Lesser General Public License for more details. 00013 00014 // You should have received a copy of the GNU Lesser General Public 00015 // License along with this library; if not, write to the Free Software 00016 // Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA 00017 00018 00019 #ifndef LIBMESH_THREADS_H 00020 #define LIBMESH_THREADS_H 00021 00022 // Local includes 00023 #include "libmesh/libmesh_config.h" 00024 #include "libmesh/libmesh_common.h" // for libmesh_assert 00025 00026 // Threading building blocks includes 00027 #ifdef LIBMESH_HAVE_TBB_API 00028 # include "libmesh/libmesh_logging.h" // only mess with the perflog if we are really multithreaded 00029 # include "tbb/tbb_stddef.h" 00030 # include "tbb/blocked_range.h" 00031 # include "tbb/parallel_for.h" 00032 # include "tbb/parallel_reduce.h" 00033 # include "tbb/task_scheduler_init.h" 00034 # include "tbb/partitioner.h" 00035 # include "tbb/spin_mutex.h" 00036 # include "tbb/recursive_mutex.h" 00037 # include "tbb/atomic.h" 00038 #endif 00039 00040 // C++ includes 00041 #ifdef LIBMESH_HAVE_STD_THREAD 00042 # include <thread> 00043 #elif LIBMESH_HAVE_TBB_CXX_THREAD 00044 # include "tbb/tbb_thread.h" 00045 #endif 00046 00047 00048 00049 namespace libMesh 00050 { 00051 00052 00057 namespace Threads 00058 { 00064 extern bool in_threads; 00065 00070 class BoolAcquire { 00071 public: 00072 explicit 00073 BoolAcquire(bool& b) : _b(b) { libmesh_assert(!_b); _b = true; } 00074 00075 ~BoolAcquire() { libmesh_assert(_b); _b = false; } 00076 private: 00077 bool& _b; 00078 }; 00079 00080 00081 #ifdef LIBMESH_HAVE_STD_THREAD 00082 //-------------------------------------------------------------------- 00086 typedef std::thread Thread; 00087 00088 #elif LIBMESH_HAVE_TBB_CXX_THREAD 00089 //-------------------------------------------------------------------- 00093 typedef tbb::tbb_thread Thread; 00094 00095 #else 00096 //-------------------------------------------------------------------- 00101 class Thread 00102 { 00103 public: 00109 template <typename Callable> 00110 Thread (Callable f) { f(); } 00111 00115 void join() {} 00116 00120 bool joinable() const { return true; } 00121 }; 00122 00123 #endif 00124 00125 00126 00127 #ifdef LIBMESH_HAVE_TBB_API 00128 //------------------------------------------------------------------- 00132 typedef tbb::task_scheduler_init task_scheduler_init; 00133 00134 00135 00136 //------------------------------------------------------------------- 00141 typedef tbb::split split; 00142 00143 00144 00145 //------------------------------------------------------------------- 00150 template <typename Range, typename Body> 00151 inline 00152 void parallel_for (const Range &range, const Body &body) 00153 { 00154 BoolAcquire b(in_threads); 00155 00156 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING 00157 const bool logging_was_enabled = libMesh::perflog.logging_enabled(); 00158 00159 if (libMesh::n_threads() > 1) 00160 libMesh::perflog.disable_logging(); 00161 #endif 00162 00163 if (libMesh::n_threads() > 1) 00164 tbb::parallel_for (range, body, tbb::auto_partitioner()); 00165 00166 else 00167 body(range); 00168 00169 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING 00170 if (libMesh::n_threads() > 1 && logging_was_enabled) 00171 libMesh::perflog.enable_logging(); 00172 #endif 00173 } 00174 00175 00176 00177 //------------------------------------------------------------------- 00182 template <typename Range, typename Body, typename Partitioner> 00183 inline 00184 void parallel_for (const Range &range, const Body &body, const Partitioner &partitioner) 00185 { 00186 BoolAcquire b(in_threads); 00187 00188 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING 00189 const bool logging_was_enabled = libMesh::perflog.logging_enabled(); 00190 00191 if (libMesh::n_threads() > 1) 00192 libMesh::perflog.disable_logging(); 00193 #endif 00194 00195 if (libMesh::n_threads() > 1) 00196 tbb::parallel_for (range, body, partitioner); 00197 00198 else 00199 body(range); 00200 00201 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING 00202 if (libMesh::n_threads() > 1 && logging_was_enabled) 00203 libMesh::perflog.enable_logging(); 00204 #endif 00205 } 00206 00207 00208 00209 //------------------------------------------------------------------- 00214 template <typename Range, typename Body> 00215 inline 00216 void parallel_reduce (const Range &range, Body &body) 00217 { 00218 BoolAcquire b(in_threads); 00219 00220 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING 00221 const bool logging_was_enabled = libMesh::perflog.logging_enabled(); 00222 00223 if (libMesh::n_threads() > 1) 00224 libMesh::perflog.disable_logging(); 00225 #endif 00226 00227 if (libMesh::n_threads() > 1) 00228 tbb::parallel_reduce (range, body, tbb::auto_partitioner()); 00229 00230 else 00231 body(range); 00232 00233 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING 00234 if (libMesh::n_threads() > 1 && logging_was_enabled) 00235 libMesh::perflog.enable_logging(); 00236 #endif 00237 } 00238 00239 00240 00241 //------------------------------------------------------------------- 00246 template <typename Range, typename Body, typename Partitioner> 00247 inline 00248 void parallel_reduce (const Range &range, Body &body, const Partitioner &partitioner) 00249 { 00250 BoolAcquire b(in_threads); 00251 00252 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING 00253 const bool logging_was_enabled = libMesh::perflog.logging_enabled(); 00254 00255 if (libMesh::n_threads() > 1) 00256 libMesh::perflog.disable_logging(); 00257 #endif 00258 00259 if (libMesh::n_threads() > 1) 00260 tbb::parallel_reduce (range, body, partitioner); 00261 00262 else 00263 body(range); 00264 00265 #ifdef LIBMESH_ENABLE_PERFORMANCE_LOGGING 00266 if (libMesh::n_threads() > 1 && logging_was_enabled) 00267 libMesh::perflog.enable_logging(); 00268 #endif 00269 } 00270 00271 00272 00273 //------------------------------------------------------------------- 00278 typedef tbb::spin_mutex spin_mutex; 00279 00280 //------------------------------------------------------------------- 00286 typedef tbb::recursive_mutex recursive_mutex; 00287 00288 //------------------------------------------------------------------- 00294 template <typename T> 00295 class atomic : public tbb::atomic<T> {}; 00296 00297 00298 00299 #else //LIBMESH_HAVE_TBB_API 00300 00301 //------------------------------------------------------------------- 00305 class task_scheduler_init 00306 { 00307 public: 00308 static const int automatic = -1; 00309 explicit task_scheduler_init (int = automatic) {} 00310 void initialize (int = automatic) {} 00311 void terminate () {} 00312 }; 00313 00314 //------------------------------------------------------------------- 00319 class split {}; 00320 00321 //------------------------------------------------------------------- 00326 template <typename Range, typename Body> 00327 inline 00328 void parallel_for (const Range &range, const Body &body) 00329 { 00330 BoolAcquire b(in_threads); 00331 body(range); 00332 } 00333 00334 //------------------------------------------------------------------- 00339 template <typename Range, typename Body, typename Partitioner> 00340 inline 00341 void parallel_for (const Range &range, const Body &body, const Partitioner &) 00342 { 00343 BoolAcquire b(in_threads); 00344 body(range); 00345 } 00346 00347 //------------------------------------------------------------------- 00352 template <typename Range, typename Body> 00353 inline 00354 void parallel_reduce (const Range &range, Body &body) 00355 { 00356 BoolAcquire b(in_threads); 00357 body(range); 00358 } 00359 00360 //------------------------------------------------------------------- 00365 template <typename Range, typename Body, typename Partitioner> 00366 inline 00367 void parallel_reduce (const Range &range, Body &body, const Partitioner &) 00368 { 00369 BoolAcquire b(in_threads); 00370 body(range); 00371 } 00372 00373 //------------------------------------------------------------------- 00378 class spin_mutex 00379 { 00380 public: 00381 spin_mutex() {} 00382 void lock () {} 00383 void unlock () {} 00384 00385 class scoped_lock 00386 { 00387 public: 00388 scoped_lock () {} 00389 explicit scoped_lock ( spin_mutex& ) {} 00390 void acquire ( spin_mutex& ) {} 00391 void release () {} 00392 }; 00393 }; 00394 00395 //------------------------------------------------------------------- 00400 class recursive_mutex 00401 { 00402 public: 00403 recursive_mutex() {} 00404 00405 class scoped_lock 00406 { 00407 public: 00408 scoped_lock () {} 00409 explicit scoped_lock ( recursive_mutex& ) {} 00410 void acquire ( recursive_mutex& ) {} 00411 void release () {} 00412 }; 00413 }; 00414 00415 //------------------------------------------------------------------- 00420 template <typename T> 00421 class atomic 00422 { 00423 public: 00424 atomic () : _val(0) {} 00425 operator T& () { return _val; } 00426 private: 00427 T _val; 00428 }; 00429 00430 00431 #endif // #ifdef LIBMESH_HAVE_TBB_API 00432 00433 00434 00438 template <typename T> 00439 class BlockedRange 00440 { 00441 public: 00445 typedef T const_iterator; 00446 00452 explicit BlockedRange (const unsigned int new_grainsize = 1000) : 00453 _grainsize(new_grainsize) 00454 {} 00455 00462 BlockedRange (const const_iterator first, 00463 const const_iterator last, 00464 const unsigned int new_grainsize = 1000) : 00465 _grainsize(new_grainsize) 00466 { 00467 this->reset(first, last); 00468 } 00469 00483 BlockedRange (const BlockedRange<T> &r): 00484 _end(r._end), 00485 _begin(r._begin), 00486 _grainsize(r._grainsize) 00487 {} 00488 00494 BlockedRange (BlockedRange<T> &r, Threads::split ) : 00495 _end(r._end), 00496 _begin(r._begin), 00497 _grainsize(r._grainsize) 00498 { 00499 const_iterator 00500 beginning = r._begin, 00501 ending = r._end, 00502 middle = beginning + (ending - beginning)/2u; 00503 00504 r._end = _begin = middle; 00505 } 00506 00510 void reset (const const_iterator first, 00511 const const_iterator last) 00512 { 00513 _begin = first; 00514 _end = last; 00515 } 00516 00520 const_iterator begin () const { return _begin; } 00521 00525 const_iterator end () const { return _end; } 00526 00531 unsigned int grainsize () const {return _grainsize;} 00532 00536 void grainsize (const unsigned int &gs) {_grainsize = gs;} 00537 00541 int size () const { return (_end -_begin); } 00542 00543 //------------------------------------------------------------------------ 00544 // Methods that implement Range concept 00545 //------------------------------------------------------------------------ 00546 00550 bool empty() const { return (_begin == _end); } 00551 00555 bool is_divisible() const { return ((_begin + this->grainsize()) < _end); } 00556 00557 private: 00558 00559 const_iterator _end; 00560 const_iterator _begin; 00561 unsigned int _grainsize; 00562 }; 00563 00564 00565 00569 extern spin_mutex spin_mtx; 00570 00574 extern recursive_mutex recursive_mtx; 00575 00576 } // namespace Threads 00577 00578 } // namespace libMesh 00579 00580 #endif // LIBMESH_THREADS_H
Site Created By: libMesh Developers
Last modified: February 05 2013 19:54:49 UTC
Hosted By: