# include <deal.II/lac/vector.h>
# include <functional>
+# include <numeric>
+
DEAL_II_NAMESPACE_OPEN
+
+namespace
+{
+ /**
+ * Optimized pack function for values assigned on degrees of freedom.
+ *
+ * Given that the elements of @p dof_values are stored in consecutive
+ * locations, we can just memcpy them. Since floating point values don't
+ * compress well, we also forgo the compression the default
+ * Utilities::pack() and Utilities::unpack() functions offer.
+ */
+ template <typename value_type>
+ std::vector<char>
+ pack_dof_values(std::vector<Vector<value_type>> &dof_values,
+ const unsigned int dofs_per_cell)
+ {
+ for (const auto &values : dof_values)
+ AssertDimension(values.size(), dofs_per_cell);
+
+ const std::size_t bytes_per_entry = sizeof(value_type) * dofs_per_cell;
+
+ std::vector<char> buffer(dof_values.size() * bytes_per_entry);
+ for (unsigned int i = 0; i < dof_values.size(); ++i)
+ std::memcpy(&buffer[i * bytes_per_entry],
+ &dof_values[i](0),
+ bytes_per_entry);
+
+ return buffer;
+ }
+
+
+
+ /**
+ * Optimized unpack function for values assigned on degrees of freedom.
+ */
+ template <typename value_type>
+ std::vector<Vector<value_type>>
+ unpack_dof_values(
+ const boost::iterator_range<std::vector<char>::const_iterator> &data_range,
+ const unsigned int dofs_per_cell)
+ {
+ const std::size_t bytes_per_entry = sizeof(value_type) * dofs_per_cell;
+ const unsigned int n_elements = data_range.size() / bytes_per_entry;
+
+ Assert((data_range.size() % bytes_per_entry == 0), ExcInternalError());
+
+ std::vector<Vector<value_type>> unpacked_data;
+ unpacked_data.reserve(n_elements);
+ for (unsigned int i = 0; i < n_elements; ++i)
+ {
+ Vector<value_type> dof_values(dofs_per_cell);
+ std::memcpy(&dof_values(0),
+ &(*std::next(data_range.begin(), i * bytes_per_entry)),
+ bytes_per_entry);
+ unpacked_data.emplace_back(std::move(dof_values));
+ }
+
+ return unpacked_data;
+ }
+} // namespace
+
+
+
namespace parallel
{
namespace distributed
typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
// create buffer for each individual object
- std::vector<::dealii::Vector<typename VectorType::value_type>> dofvalues(
+ std::vector<::dealii::Vector<typename VectorType::value_type>> dof_values(
input_vectors.size());
+ unsigned int fe_index = 0;
if (DoFHandlerType::is_hp_dof_handler)
{
- unsigned int fe_index = numbers::invalid_unsigned_int;
switch (status)
{
case parallel::distributed::Triangulation<
Assert(false, ExcInternalError());
break;
}
+ }
- const unsigned int dofs_per_cell =
- dof_handler->get_fe(fe_index).dofs_per_cell;
+ const unsigned int dofs_per_cell =
+ dof_handler->get_fe(fe_index).dofs_per_cell;
- auto it_input = input_vectors.cbegin();
- auto it_output = dofvalues.begin();
- for (; it_input != input_vectors.cend(); ++it_input, ++it_output)
- {
- it_output->reinit(dofs_per_cell);
- cell->get_interpolated_dof_values(*(*it_input),
- *it_output,
- fe_index);
- }
- }
- else
+ auto it_input = input_vectors.cbegin();
+ auto it_output = dof_values.begin();
+ for (; it_input != input_vectors.cend(); ++it_input, ++it_output)
{
- auto it_input = input_vectors.cbegin();
- auto it_output = dofvalues.begin();
- for (; it_input != input_vectors.cend(); ++it_input, ++it_output)
- {
- it_output->reinit(cell->get_fe().dofs_per_cell);
- cell->get_interpolated_dof_values(*(*it_input), *it_output);
- }
+ it_output->reinit(dofs_per_cell);
+ cell->get_interpolated_dof_values(*(*it_input), *it_output, fe_index);
}
- // We choose to compress the packed data depending on the registered
- // dofhandler object. In case of hp, we write in the variable size buffer
- // and thus allow compression. In the other case, we use the fixed size
- // buffer and require consistent data sizes on each cell, so we leave it
- // uncompressed.
- return Utilities::pack(
- dofvalues, /*allow_compression=*/DoFHandlerType::is_hp_dof_handler);
+ return pack_dof_values<typename VectorType::value_type>(dof_values,
+ dofs_per_cell);
}
{
typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
- const std::vector<::dealii::Vector<typename VectorType::value_type>>
- dofvalues = Utilities::unpack<
- std::vector<::dealii::Vector<typename VectorType::value_type>>>(
- data_range.begin(),
- data_range.end(),
- /*allow_compression=*/DoFHandlerType::is_hp_dof_handler);
-
- // check if sizes match
- Assert(dofvalues.size() == all_out.size(), ExcInternalError());
-
+ unsigned int fe_index = 0;
if (DoFHandlerType::is_hp_dof_handler)
{
- unsigned int fe_index = numbers::invalid_unsigned_int;
-
switch (status)
{
case parallel::distributed::Triangulation<
Assert(false, ExcInternalError());
break;
}
-
- // check if we have enough dofs provided by the FE object
- // to interpolate the transferred data correctly
- for (auto it_dofvalues = dofvalues.begin();
- it_dofvalues != dofvalues.end();
- ++it_dofvalues)
- Assert(
- dof_handler->get_fe(fe_index).dofs_per_cell ==
- it_dofvalues->size(),
- ExcMessage(
- "The transferred data was packed with a different number of dofs than the "
- "currently registered FE object assigned to the DoFHandler has."));
-
- // distribute data for each registered vector on mesh
- auto it_input = dofvalues.cbegin();
- auto it_output = all_out.begin();
- for (; it_input != dofvalues.cend(); ++it_input, ++it_output)
- cell->set_dof_values_by_interpolation(*it_input,
- *(*it_output),
- fe_index);
}
- else
- {
- // check if we have enough dofs provided by the FE object
- // to interpolate the transferred data correctly
- for (auto it_dofvalues = dofvalues.begin();
- it_dofvalues != dofvalues.end();
- ++it_dofvalues)
- Assert(
- cell->get_fe().dofs_per_cell == it_dofvalues->size(),
- ExcMessage(
- "The transferred data was packed with a different number of dofs than the "
- "currently registered FE object assigned to the DoFHandler has."));
-
- // distribute data for each registered vector on mesh
- auto it_input = dofvalues.cbegin();
- auto it_output = all_out.begin();
- for (; it_input != dofvalues.cend(); ++it_input, ++it_output)
- cell->set_dof_values_by_interpolation(*it_input, *(*it_output));
- }
- }
+ const unsigned int dofs_per_cell =
+ dof_handler->get_fe(fe_index).dofs_per_cell;
+ const std::vector<::dealii::Vector<typename VectorType::value_type>>
+ dof_values =
+ unpack_dof_values<typename VectorType::value_type>(data_range,
+ dofs_per_cell);
+
+ // check if sizes match
+ Assert(dof_values.size() == all_out.size(), ExcInternalError());
+
+ // check if we have enough dofs provided by the FE object
+ // to interpolate the transferred data correctly
+ for (auto it_dof_values = dof_values.begin();
+ it_dof_values != dof_values.end();
+ ++it_dof_values)
+ Assert(
+ dofs_per_cell == it_dof_values->size(),
+ ExcMessage(
+ "The transferred data was packed with a different number of dofs than the "
+ "currently registered FE object assigned to the DoFHandler has."));
+
+ // distribute data for each registered vector on mesh
+ auto it_input = dof_values.cbegin();
+ auto it_output = all_out.begin();
+ for (; it_input != dof_values.cend(); ++it_input, ++it_output)
+ cell->set_dof_values_by_interpolation(*it_input,
+ *(*it_output),
+ fe_index);
+ }
} // namespace distributed
} // namespace parallel