<h3>Specific improvements</h3>
<ol>
+ <li> Improved: MPI collective operations such as MPI::sum, MPI::max now
+ check for job_supports_mpi() internally, which allows running them also
+ without a call to MPI_Init.
+ <br>
+ (Martin Kronbichler, 2015/01/13)
+ </li>
+
+ <li> Changed: The method job_supports_mpi() now resides in the namespace
+ Utilities::MPI instead of Utilities::System for consistency with other MPI
+ methods. The old method has been marked deprecated and will be removed in
+ a future version.
+ <br>
+ (Martin Kronbichler, 2015/01/13)
+ </li>
+
<li> Fixed: The update of ghost values in parallel::distributed::Vector when
calling the assignment operator is now active when one of the two vector had
its ghost values updated before or when the layout of the right hand side
// ---------------------------------------------------------------------
//
-// Copyright (C) 2011 - 2014 by the deal.II authors
+// Copyright (C) 2011 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
~MPI_InitFinalize();
};
+ /**
+ * Return whether (i) deal.II has been compiled to support MPI (for
+ * example by compiling with <code>CXX=mpiCC</code>) and if so whether
+ * (ii) <code>MPI_Init()</code> has been called (for example using the
+ * Utilities::MPI::MPI_InitFinalize class). In other words, the result
+ * indicates whether the current job is running under MPI.
+ *
+ * @note The function does not take into account whether an MPI job
+ * actually runs on more than one processor or is, in fact, a single-node
+ * job that happens to run under MPI.
+ */
+ bool job_supports_mpi ();
+
namespace internal
{
#ifdef DEAL_II_WITH_MPI
const MPI_Comm &mpi_communicator)
{
#ifdef DEAL_II_WITH_MPI
- T sum;
- MPI_Allreduce (const_cast<void *>(static_cast<const void *>(&t)),
- &sum, 1, internal::mpi_type_id(&t), MPI_SUM,
- mpi_communicator);
- return sum;
-#else
- (void)mpi_communicator;
- return t;
+ if (job_supports_mpi())
+ {
+ T sum;
+ MPI_Allreduce (const_cast<void *>(static_cast<const void *>(&t)),
+ &sum, 1, internal::mpi_type_id(&t), MPI_SUM,
+ mpi_communicator);
+ return sum;
+ }
+ else
#endif
+ {
+ (void)mpi_communicator;
+ return t;
+ }
}
T (&sums)[N])
{
#ifdef DEAL_II_WITH_MPI
- MPI_Allreduce ((&values[0] != &sums[0]
- ?
- const_cast<void *>(static_cast<const void *>(&values[0]))
- :
- MPI_IN_PLACE),
- &sums[0], N, internal::mpi_type_id(values), MPI_SUM,
- mpi_communicator);
-#else
- (void)mpi_communicator;
- for (unsigned int i=0; i<N; ++i)
- sums[i] = values[i];
+ if (job_supports_mpi())
+ {
+ MPI_Allreduce ((&values[0] != &sums[0]
+ ?
+ const_cast<void *>(static_cast<const void *>(&values[0]))
+ :
+ MPI_IN_PLACE),
+ &sums[0], N, internal::mpi_type_id(values), MPI_SUM,
+ mpi_communicator);
+ }
+ else
#endif
+ {
+ (void)mpi_communicator;
+ for (unsigned int i=0; i<N; ++i)
+ sums[i] = values[i];
+ }
}
std::vector<T> &sums)
{
#ifdef DEAL_II_WITH_MPI
- sums.resize (values.size());
- MPI_Allreduce ((&values[0] != &sums[0]
- ?
- const_cast<void *>(static_cast<const void *>(&values[0]))
- :
- MPI_IN_PLACE),
- &sums[0], values.size(), internal::mpi_type_id((T *)0), MPI_SUM,
- mpi_communicator);
-#else
- (void)mpi_communicator;
- sums = values;
+ if (job_supports_mpi())
+ {
+ sums.resize (values.size());
+ MPI_Allreduce ((&values[0] != &sums[0]
+ ?
+ const_cast<void *>(static_cast<const void *>(&values[0]))
+ :
+ MPI_IN_PLACE),
+ &sums[0], values.size(), internal::mpi_type_id((T *)0), MPI_SUM,
+ mpi_communicator);
+ }
+ else
#endif
+ {
+ (void)mpi_communicator;
+ sums = values;
+ }
}
const MPI_Comm &mpi_communicator)
{
#ifdef DEAL_II_WITH_MPI
- T sum;
- MPI_Allreduce (const_cast<void *>(static_cast<const void *>(&t)),
- &sum, 1, internal::mpi_type_id(&t), MPI_MAX,
- mpi_communicator);
- return sum;
-#else
- (void)mpi_communicator;
- return t;
+ if (job_supports_mpi())
+ {
+ T sum;
+ MPI_Allreduce (const_cast<void *>(static_cast<const void *>(&t)),
+ &sum, 1, internal::mpi_type_id(&t), MPI_MAX,
+ mpi_communicator);
+ return sum;
+ }
+ else
#endif
+ {
+ (void)mpi_communicator;
+ return t;
+ }
}
T (&maxima)[N])
{
#ifdef DEAL_II_WITH_MPI
- MPI_Allreduce ((&values[0] != &maxima[0]
- ?
- const_cast<void *>(static_cast<const void *>(&values[0]))
- :
- MPI_IN_PLACE),
- &maxima[0], N, internal::mpi_type_id(values), MPI_MAX,
- mpi_communicator);
-#else
- (void)mpi_communicator;
- for (unsigned int i=0; i<N; ++i)
- maxima[i] = values[i];
+ if (job_supports_mpi())
+ {
+ MPI_Allreduce ((&values[0] != &maxima[0]
+ ?
+ const_cast<void *>(static_cast<const void *>(&values[0]))
+ :
+ MPI_IN_PLACE),
+ &maxima[0], N, internal::mpi_type_id(values), MPI_MAX,
+ mpi_communicator);
+ }
+ else
#endif
+ {
+ (void)mpi_communicator;
+ for (unsigned int i=0; i<N; ++i)
+ maxima[i] = values[i];
+ }
}
std::vector<T> &maxima)
{
#ifdef DEAL_II_WITH_MPI
- maxima.resize (values.size());
- MPI_Allreduce ((&values[0] != &maxima[0]
- ?
- const_cast<void *>(static_cast<const void *>(&values[0]))
- :
- MPI_IN_PLACE),
- &maxima[0], values.size(), internal::mpi_type_id((T *)0), MPI_MAX,
- mpi_communicator);
+ if (job_supports_mpi())
+ {
+ maxima.resize (values.size());
+ MPI_Allreduce ((&values[0] != &maxima[0]
+ ?
+ const_cast<void *>(static_cast<const void *>(&values[0]))
+ :
+ MPI_IN_PLACE),
+ &maxima[0], values.size(), internal::mpi_type_id((T *)0), MPI_MAX,
+ mpi_communicator);
+ }
+ else
+#endif
+ {
+ (void)mpi_communicator;
+ maxima = values;
+ }
+ }
+
+
+ inline
+ bool job_supports_mpi ()
+ {
+#ifdef DEAL_II_WITH_MPI
+ int MPI_has_been_started = 0;
+ MPI_Initialized(&MPI_has_been_started);
+
+ return (MPI_has_been_started > 0);
#else
- (void)mpi_communicator;
- maxima = values;
+ return false;
#endif
}
} // end of namespace MPI
// ---------------------------------------------------------------------
//
-// Copyright (C) 2005 - 2013 by the deal.II authors
+// Copyright (C) 2005 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
std::string get_time ();
/**
- * Return whether (i) deal.II has been compiled to support MPI (for
- * example by compiling with <code>CXX=mpiCC</code>) and if so whether
- * (ii) <code>MPI_Init()</code> has been called (for example using the
- * Utilities::System::MPI_InitFinalize class). In other words, the result
- * indicates whether the current job is running under MPI.
- *
- * @note The function does not take into account whether an MPI job
- * actually runs on more than one processor or is, in fact, a single-node
- * job that happens to run under MPI.
+ * @deprecated Use Utilities::MPI::job_supports_mpi() instead.
*/
- bool job_supports_mpi ();
+ bool job_supports_mpi () DEAL_II_DEPRECATED;
}
#ifdef DEAL_II_WITH_MPI
#ifdef DEBUG
- if (Utilities::System::job_supports_mpi())
+ if (Utilities::MPI::job_supports_mpi())
{
// make sure that there are not outstanding requests from updating
// ghost values or compress
// ---------------------------------------------------------------------
//
-// Copyright (C) 2011 - 2014 by the deal.II authors
+// Copyright (C) 2011 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
dynamic_cast<const parallel::distributed::Triangulation<dim>*>(&tria);
if (dist_tria != 0)
{
- if (Utilities::System::job_supports_mpi())
+ if (Utilities::MPI::job_supports_mpi())
{
int communicators_same = 0;
MPI_Comm_compare (dist_tria->get_communicator(), comm_mf,
internal::assert_communicator_equality (dof_handler[0]->get_tria(),
additional_data.mpi_communicator);
size_info.communicator = additional_data.mpi_communicator;
- if (Utilities::System::job_supports_mpi() == true)
+ if (Utilities::MPI::job_supports_mpi() == true)
{
size_info.my_pid =
Utilities::MPI::this_mpi_process(size_info.communicator);
internal::assert_communicator_equality (dof_handler[0]->get_tria(),
additional_data.mpi_communicator);
size_info.communicator = additional_data.mpi_communicator;
- if (Utilities::System::job_supports_mpi() == true)
+ if (Utilities::MPI::job_supports_mpi() == true)
{
size_info.my_pid =
Utilities::MPI::this_mpi_process(size_info.communicator);
void SizeInfo::print_memory_statistics (STREAM &out,
std::size_t data_length) const
{
- Utilities::MPI::MinMaxAvg memory_c;
- if (Utilities::System::job_supports_mpi() == true)
- {
- memory_c = Utilities::MPI::min_max_avg (1e-6*data_length,
- communicator);
- }
- else
- {
- memory_c.sum = 1e-6*data_length;
- memory_c.min = memory_c.sum;
- memory_c.max = memory_c.sum;
- memory_c.avg = memory_c.sum;
- memory_c.min_index = 0;
- memory_c.max_index = 0;
- }
+ Utilities::MPI::MinMaxAvg memory_c
+ = Utilities::MPI::min_max_avg (1e-6*data_length, communicator);
if (n_procs < 2)
out << memory_c.min;
else
// ---------------------------------------------------------------------
//
-// Copyright (C) 2005 - 2014 by the deal.II authors
+// Copyright (C) 2005 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
min_max_avg(const double my_value,
const MPI_Comm &mpi_communicator)
{
+ // If MPI was not started, we have a serial computation and cannot run
+ // the other MPI commands
+ if (job_supports_mpi() == false)
+ {
+ MinMaxAvg result;
+ result.sum = my_value;
+ result.avg = my_value;
+ result.min = my_value;
+ result.max = my_value;
+ result.min_index = 0;
+ result.max_index = 0;
+
+ return result;
+ }
+
// To avoid uninitialized values on some MPI implementations, provide
// result with a default value already...
MinMaxAvg result = { 0., std::numeric_limits<double>::max(),
// when running PETSc, because we initialize MPI ourselves before calling
// PetscInitialize
#ifdef DEAL_II_WITH_MPI
- int MPI_has_been_started = 0;
- MPI_Initialized(&MPI_has_been_started);
- if (Utilities::System::job_supports_mpi() == true &&
- MPI_has_been_started != 0)
+ if (job_supports_mpi() == true)
{
if (std::uncaught_exception())
{
n_ghost_indices_data = ghost_indices_data.n_elements();
have_ghost_indices =
- Utilities::System::job_supports_mpi() &&
Utilities::MPI::sum(n_ghost_indices_data, communicator) > 0;
// In the rest of this function, we determine the point-to-point
// ---------------------------------------------------------------------
//
-// Copyright (C) 1998 - 2014 by the deal.II authors
+// Copyright (C) 1998 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
#endif
#ifdef DEAL_II_WITH_MPI
- if (sync_wall_time && Utilities::System::job_supports_mpi())
+ if (sync_wall_time && Utilities::MPI::job_supports_mpi())
{
this->mpi_data
= Utilities::MPI::min_max_avg (last_lap_time, mpi_communicator);
double Timer::get_lap_time() const
{
- // time already has the difference
- // between the last start()/stop()
- // cycle.
-#ifdef DEAL_II_WITH_MPI
- if (Utilities::System::job_supports_mpi())
- return Utilities::MPI::max (last_lap_time, mpi_communicator);
- else
-#endif
- return last_lap_time;
-
+ // time already has the difference between the last start()/stop() cycle.
+ return Utilities::MPI::max (last_lap_time, mpi_communicator);
}
const double running_time = dtime - start_time + dtime_children
- start_time_children + cumulative_time;
- if (Utilities::System::job_supports_mpi())
- // in case of MPI, need to get the time
- // passed by summing the time over all
- // processes in the network. works also
- // in case we just want to have the time
- // of a single thread, since then the
- // communicator is MPI_COMM_SELF
- return Utilities::MPI::sum (running_time, mpi_communicator);
- else
- return running_time;
+ // in case of MPI, need to get the time passed by summing the time over
+ // all processes in the network. works also in case we just want to have
+ // the time of a single thread, since then the communicator is
+ // MPI_COMM_SELF
+ return Utilities::MPI::sum (running_time, mpi_communicator);
+
#elif defined(DEAL_II_MSVC)
const double running_time = windows::cpu_clock() - start_time + cumulative_time;
return running_time;
}
else
{
- if (Utilities::System::job_supports_mpi())
- return Utilities::MPI::sum (cumulative_time, mpi_communicator);
- else
- return cumulative_time;
+ return Utilities::MPI::sum (cumulative_time, mpi_communicator);
}
}
// initialize the Timers here according to that
double cpu_time = sections[actual_section_name].timer();
sections[actual_section_name].total_cpu_time
- += (Utilities::System::job_supports_mpi()
- ?
- Utilities::MPI::sum (cpu_time, mpi_communicator)
- :
- cpu_time);
+ += Utilities::MPI::sum (cpu_time, mpi_communicator);
// in case we have to print out something, do that here...
if ((output_frequency == every_call || output_frequency == every_call_and_summary)
// in case we want to write CPU times
if (output_type != wall_times)
{
- double total_cpu_time = (Utilities::System::job_supports_mpi()
- ?
- Utilities::MPI::sum (timer_all(),
- mpi_communicator)
- :
- timer_all());
+ double total_cpu_time = Utilities::MPI::sum(timer_all(), mpi_communicator);
// check that the sum of all times is
// less or equal than the total
// ---------------------------------------------------------------------
//
-// Copyright (C) 2005 - 2013 by the deal.II authors
+// Copyright (C) 2005 - 2015 by the deal.II authors
//
// This file is part of the deal.II library.
//
bool job_supports_mpi ()
{
-#ifdef DEAL_II_WITH_MPI
- int MPI_has_been_started = 0;
- MPI_Initialized(&MPI_has_been_started);
-
- return true && (MPI_has_been_started > 0);
-#else
- return false;
-#endif
+ return Utilities::MPI::job_supports_mpi();
}
}
void test()
{
- Assert( Utilities::System::job_supports_mpi(), ExcInternalError());
+ Assert( Utilities::MPI::job_supports_mpi(), ExcInternalError());
unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
void test_mpi()
{
- Assert( Utilities::System::job_supports_mpi(), ExcInternalError());
+ Assert( Utilities::MPI::job_supports_mpi(), ExcInternalError());
unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
void test_mpi()
{
- Assert( Utilities::System::job_supports_mpi(), ExcInternalError());
+ Assert( Utilities::MPI::job_supports_mpi(), ExcInternalError());
unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
void test_mpi()
{
- Assert( Utilities::System::job_supports_mpi(), ExcInternalError());
+ Assert( Utilities::MPI::job_supports_mpi(), ExcInternalError());
unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);
const unsigned int numprocs = Utilities::MPI::n_mpi_processes (MPI_COMM_WORLD);
void test_mpi()
{
- Assert( Utilities::System::job_supports_mpi(), ExcInternalError());
+ Assert( Utilities::MPI::job_supports_mpi(), ExcInternalError());
unsigned int myid = Utilities::MPI::this_mpi_process (MPI_COMM_WORLD);