<ol>
+ <li> New: MultithreadInfo::set_thread_limit() can now be called more than
+ once and the environment variable DEAL_II_NUM_THREADS will be respected
+ even if user code never calls it.
+ <br>
+ (Timo Heister, 2015/07/26)
+ </li>
+
<li> New: IndexSet now implements iterators.
<br>
(Timo Heister, 2015/07/12)
* given, the default from TBB is used (based on the number of cores in the
* system).
*
- * Due to limitations in the way TBB can be controlled, only the first call
- * to this method will have any effect. In practice, this means that you
- * need to call this function before you get to any point in your program
- * where multiple threads may be created. In other words, the correct place
- * for a call to this function is at the top of your <code>main()</code>
- * function.
- *
- * This routine is called automatically by MPI_InitFinalize. Use the
- * appropriate argument of the constructor of MPI_InitFinalize if you have
- * an MPI based code.
+ * This routine is executed automatically with the default argument before
+ * your code in main() is running (using a static constructor). It is also
+ * executed by MPI_InitFinalize. Use the appropriate argument of the
+ * constructor of MPI_InitFinalize if you have an MPI based code.
*/
static void set_thread_limit (const unsigned int max_threads = numbers::invalid_unsigned_int);
-
/**
* Returns if the TBB is running using a single thread either because of
* thread affinity or because it is set via a call to set_thread_limit. This
void MultithreadInfo::set_thread_limit(const unsigned int max_threads)
{
- Assert(n_max_threads==numbers::invalid_unsigned_int,
- ExcMessage("Calling set_thread_limit() more than once is not supported!"));
-
// set the maximal number of threads to the given value as specified
n_max_threads = max_threads;
n_max_threads = max_threads_env;
}
}
- // finally see if we need to tell the TBB about this. if no value has been set
- // (the if-branch), then simply set n_max_threads to what we get from the TBB
- // but don't do anything further. otherwise (some value has been set), let the
- // TBB use this value
+ // Without restrictions from the user query TBB for the recommended number
+ // of threads:
if (n_max_threads == numbers::invalid_unsigned_int)
n_max_threads = tbb::task_scheduler_init::default_num_threads();
- else
- {
- static tbb::task_scheduler_init dummy (n_max_threads);
- }
+
+ // Initialize the scheduler and destroy the old one before doing so
+ static tbb::task_scheduler_init dummy (tbb::task_scheduler_init::deferred);
+ if (dummy.is_active())
+ dummy.terminate();
+ dummy.initialize(n_max_threads);
}
unsigned int MultithreadInfo::n_threads()
{
- if (n_max_threads == numbers::invalid_unsigned_int)
- return tbb::task_scheduler_init::default_num_threads();
- else
- return n_max_threads;
+ Assert(n_max_threads != numbers::invalid_unsigned_int, ExcInternalError());
+ return n_max_threads;
}
// definition of the variable which is declared `extern' in the .h file
MultithreadInfo multithread_info;
+namespace
+{
+// Force the first call to set_thread_limit happen before any tasks in TBB are
+// used. This is necessary as tbb::task_scheduler_init has no effect if TBB
+// got automatically initialized (which happens the first time we use it).
+ struct DoOnce
+ {
+ DoOnce ()
+ {
+ MultithreadInfo::set_thread_limit (numbers::invalid_unsigned_int);
+ }
+ } do_once;
+}
DEAL_II_NAMESPACE_CLOSE
--- /dev/null
+// ---------------------------------------------------------------------
+//
+// Copyright (C) 2009 - 2015 by the deal.II authors
+//
+// This file is part of the deal.II library.
+//
+// The deal.II library is free software; you can use it, redistribute
+// it, and/or modify it under the terms of the GNU Lesser General
+// Public License as published by the Free Software Foundation; either
+// version 2.1 of the License, or (at your option) any later version.
+// The full text of the license can be found in the file LICENSE at
+// the top level of the deal.II distribution.
+//
+// ---------------------------------------------------------------------
+
+// test calling MultithreadInfo::set_thread_limit() multiple times and the
+// effect for running tasks
+
+#include "../tests.h"
+#include <iomanip>
+#include <fstream>
+#include <unistd.h>
+
+#include <deal.II/base/thread_management.h>
+
+Threads::Mutex mutex;
+unsigned int n_running;
+unsigned int n_max_running;
+
+void a_task ()
+{
+ mutex.acquire();
+ n_running++;
+ n_max_running=std::max(n_max_running, n_running);
+ mutex.release();
+
+ // Sleep some time to make sure all other concurrently running tasks enter
+ // here.
+ sleep (1);
+
+ mutex.acquire();
+ n_running--;
+ mutex.release();
+}
+
+
+void test()
+{
+ Threads::TaskGroup<> tg;
+
+ n_running = 0;
+ n_max_running = 0;
+
+ // force all tasks to wait until we are done starting
+ mutex.acquire();
+
+ for (unsigned int t=0; t<10; ++t)
+ tg += Threads::new_task(a_task);
+
+ // now let the tasks run
+ mutex.release();
+
+ tg.join_all();
+ deallog << "max concurrent running: " << n_max_running
+ << " should be: " << MultithreadInfo::n_threads()
+ << std::endl;
+}
+
+
+int main()
+{
+ std::ofstream logfile("output");
+ deallog.attach(logfile);
+ deallog.depth_console(0);
+ deallog.threshold_double(1.e-10);
+
+ const unsigned int n = testing_max_num_threads();
+
+ // if we have a machine with less than 5 cores fake the output:
+ if (MultithreadInfo::n_threads()<n)
+ {
+ // you should buy more cores, seriously.
+ deallog << "* default thread limit for tests: " << n << std::endl;
+ deallog << "max concurrent running: " << n
+ << " should be: " << n
+ << std::endl;
+ }
+ else if (MultithreadInfo::n_threads()==n)
+ {
+ deallog << "* default thread limit for tests: " << MultithreadInfo::n_threads() << std::endl;
+ test();
+ }
+ else
+ Assert(false, ExcInternalError("did somebody mess with LimitConcurrency in tests.h?"));
+
+ deallog << "* now with thread limit 1:" << std::endl;
+ MultithreadInfo::set_thread_limit(1);
+ test();
+ deallog << "* now with thread limit 2:" << std::endl;
+ MultithreadInfo::set_thread_limit(2);
+ test();
+
+ deallog << "OK" << std::endl;
+}
--- /dev/null
+
+DEAL::* default thread limit for tests: 5
+DEAL::max concurrent running: 5 should be: 5
+DEAL::* now with thread limit 1:
+DEAL::max concurrent running: 1 should be: 1
+DEAL::* now with thread limit 2:
+DEAL::max concurrent running: 2 should be: 2
+DEAL::OK
// in real space in two locations inside the cell.
// To make sure the cell similarity code is used, only run the program with
-// one thread. For this we need to disable LimitConcurreny in ../tests.h and
-// then set the number of threads to 1 (in main() below):
-#define DEAL_II_TEST_DO_NOT_SET_THREAD_LIMIT
+// one thread.
#include "../tests.h"
#include <deal.II/base/logstream.h>
* each of them runs 64 threads, then we get astronomical loads.
* Limit concurrency to a fixed (small) number of threads, independent
* of the core count.
- *
- * Note that we can't do this if we run in MPI mode because then
- * MPI_InitFinalize already calls this function. Since every test
- * calls MPI_InitFinalize itself, we can't adjust the thread count
- * for this here.
*/
-#if !defined(DEAL_II_WITH_MPI) && !defined(DEAL_II_TEST_DO_NOT_SET_THREAD_LIMIT)
+inline unsigned int testing_max_num_threads()
+{
+ return 5;
+}
+
struct LimitConcurrency
{
LimitConcurrency ()
{
- MultithreadInfo::set_thread_limit (5);
+ MultithreadInfo::set_thread_limit (testing_max_num_threads());
}
} limit_concurrency;
-#endif