Pack functions return the packed data.
* objects that can later be retrieved after refinement, coarsening and
* repartitioning.
*/
- void
+ std::vector<char>
pack_function(
const typename parallel::distributed::Triangulation<dim>::cell_iterator
&cell,
const typename parallel::distributed::Triangulation<dim>::CellStatus
- status,
- std::vector<char> &data);
+ status);
/**
* A callback function used to unpack the data on the current mesh that
const typename parallel::distributed::Triangulation<dim>::CellStatus
status,
const boost::iterator_range<std::vector<char>::const_iterator>
- data_range);
+ &data_range);
/**
* FiniteElement used to project data from and to quadrature points.
&ContinuousQuadratureDataTransfer<dim, DataType>::pack_function,
this,
std::placeholders::_1,
- std::placeholders::_2,
- std::placeholders::_3));
+ std::placeholders::_2));
}
template <int dim, typename DataType>
- void
+ std::vector<char>
ContinuousQuadratureDataTransfer<dim, DataType>::pack_function(
const typename parallel::distributed::Triangulation<dim>::cell_iterator
&cell,
const typename parallel::distributed::Triangulation<
- dim>::CellStatus /*status*/,
- std::vector<char> &data)
+ dim>::CellStatus /*status*/)
{
pack_cell_data(cell, data_storage, matrix_quadrature);
// to get consistent data sizes on each cell for the fixed size transfer,
// we won't allow compression
- Utilities::pack(matrix_dofs, data, /*allow_compression=*/false);
+ return Utilities::pack(matrix_dofs, /*allow_compression=*/false);
}
const typename parallel::distributed::Triangulation<dim>::cell_iterator
&cell,
const typename parallel::distributed::Triangulation<dim>::CellStatus
- status,
- const boost::iterator_range<std::vector<char>::const_iterator> data_range)
+ status,
+ const boost::iterator_range<std::vector<char>::const_iterator>
+ &data_range)
{
Assert((status !=
parallel::distributed::Triangulation<dim, dim>::CELL_COARSEN),
* objects that can later be retrieved after refinement, coarsening and
* repartitioning.
*/
- void
+ std::vector<char>
pack_callback(
const typename Triangulation<dim, DoFHandlerType::space_dimension>::
cell_iterator &cell,
const typename Triangulation<dim, DoFHandlerType::space_dimension>::
- CellStatus status,
- std::vector<char> &data);
+ CellStatus status);
/**
* A callback function used to unpack the data on the current mesh that
const typename Triangulation<dim, DoFHandlerType::space_dimension>::
CellStatus status,
const boost::iterator_range<std::vector<char>::const_iterator>
- data_range,
+ & data_range,
std::vector<VectorType *> &all_out);
* of classes that do this is parallel::distributed::SolutionTransfer
* where each parallel::distributed::SolutionTransfer object that works
* on the current Triangulation object then needs to register its intent.
- * Each of these parties registers a call-back function (the second
+ * Each of these parties registers a callback function (the second
* argument here, @p pack_callback) that will be called whenever the
* triangulation's execute_coarsening_and_refinement() or save()
* functions are called.
* using save() and load(), then the cell status argument with which
* the callback is called will always be `CELL_PERSIST`.
*
- * The third argument given to the callback is a reference to the
- * buffer on which the packed data will be appended at the back.
+ * The callback function is expected to return a memory chunk of the
+ * format `std::vector<char>`, representing the packed data on a
+ * certain cell.
*
* @note The purpose of this function is to register intent to
* attach data for a single, subsequent call to
*/
unsigned int
register_data_attach(
- const std::function<void(const cell_iterator &,
- const CellStatus,
- std::vector<char> &)> &pack_callback);
+ const std::function<std::vector<char>(const cell_iterator &,
+ const CellStatus)>
+ &pack_callback);
/**
* This function is the opposite of register_data_attach(). It is called
void
notify_ready_to_unpack(
const unsigned int handle,
- const std::function<
- void(const cell_iterator &,
- const CellStatus,
- const boost::iterator_range<std::vector<char>::const_iterator>)>
+ const std::function<void(
+ const cell_iterator &,
+ const CellStatus,
+ const boost::iterator_range<std::vector<char>::const_iterator> &)>
&unpack_callback);
/**
*/
unsigned int n_attached_deserialize;
- using pack_callback_t = std::function<
- void(typename Triangulation<dim, spacedim>::cell_iterator,
- CellStatus,
- std::vector<char> &)>;
+ using pack_callback_t = std::function<std::vector<char>(
+ typename Triangulation<dim, spacedim>::cell_iterator,
+ CellStatus)>;
/**
* These callback functions will be stored in the order on how they have
const typename dealii::Triangulation<dim, spacedim>::cell_iterator
&,
const typename dealii::Triangulation<dim, spacedim>::CellStatus &,
- const boost::iterator_range<std::vector<char>::const_iterator>)>
+ const boost::iterator_range<std::vector<char>::const_iterator> &)>
&unpack_callback) const;
/**
*/
unsigned int
register_data_attach(
- const std::function<void(
+ const std::function<std::vector<char>(
const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
- const typename dealii::Triangulation<1, spacedim>::CellStatus,
- std::vector<char> &)> &pack_callback);
+ const typename dealii::Triangulation<1, spacedim>::CellStatus)>
+ &pack_callback);
/**
* This function is not implemented, but needs to be present for the
const std::function<void(
const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
const typename dealii::Triangulation<1, spacedim>::CellStatus,
- const boost::iterator_range<std::vector<char>::const_iterator>)>
+ const boost::iterator_range<std::vector<char>::const_iterator> &)>
&unpack_callback);
/**
* before a refinement step. All particles have to be attached to their
* cell to be sent around to the new processes.
*/
- void
+ std::vector<char>
store_particles(
const typename Triangulation<dim, spacedim>::cell_iterator &cell,
- const typename Triangulation<dim, spacedim>::CellStatus status,
- std::vector<char> & data) const;
+ const typename Triangulation<dim, spacedim>::CellStatus status) const;
/**
* Called by listener functions after a refinement step. The local map
const typename Triangulation<dim, spacedim>::cell_iterator &cell,
const typename Triangulation<dim, spacedim>::CellStatus status,
const boost::iterator_range<std::vector<char>::const_iterator>
- data_range);
+ &data_range);
};
/* ---------------------- inline and template functions ------------------ */
&SolutionTransfer<dim, VectorType, DoFHandlerType>::pack_callback,
this,
std::placeholders::_1,
- std::placeholders::_2,
- std::placeholders::_3));
+ std::placeholders::_2));
}
template <int dim, typename VectorType, typename DoFHandlerType>
- void
+ std::vector<char>
SolutionTransfer<dim, VectorType, DoFHandlerType>::pack_callback(
const typename Triangulation<dim, DoFHandlerType::space_dimension>::
cell_iterator &cell_,
const typename Triangulation<dim, DoFHandlerType::space_dimension>::
- CellStatus /*status*/,
- std::vector<char> &data)
+ CellStatus /*status*/)
{
typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
// to get consistent data sizes on each cell for the fixed size transfer,
// we won't allow compression
- Utilities::pack(dofvalues, data, /*allow_compression=*/false);
+ return Utilities::pack(dofvalues, /*allow_compression=*/false);
}
cell_iterator &cell_,
const typename Triangulation<dim, DoFHandlerType::space_dimension>::
CellStatus /*status*/,
- const boost::iterator_range<std::vector<char>::const_iterator> data_range,
- std::vector<VectorType *> & all_out)
+ const boost::iterator_range<std::vector<char>::const_iterator>
+ & data_range,
+ std::vector<VectorType *> &all_out)
{
typename DoFHandlerType::cell_iterator cell(*cell_, dof_handler);
callback_it != pack_callbacks.cend();
++callback_it, ++data_it)
{
- (*callback_it)(dealii_cell, cell_status, *data_it);
+ *data_it = (*callback_it)(dealii_cell, cell_status);
}
}
} // loop over quad_cell_relations
const std::function<void(
const typename dealii::Triangulation<dim, spacedim>::cell_iterator &,
const typename dealii::Triangulation<dim, spacedim>::CellStatus &,
- const boost::iterator_range<std::vector<char>::const_iterator>)>
+ const boost::iterator_range<std::vector<char>::const_iterator> &)>
&unpack_callback) const
{
Assert(cumulative_sizes_fixed.size() > 0,
template <int dim, int spacedim>
unsigned int
Triangulation<dim, spacedim>::register_data_attach(
- const std::function<void(const cell_iterator &,
- const CellStatus,
- std::vector<char> &)> &pack_callback)
+ const std::function<std::vector<char>(const cell_iterator &,
+ const CellStatus)> &pack_callback)
{
// add new callback_set
cell_attached_data.pack_callbacks.push_back(pack_callback);
const std::function<
void(const cell_iterator &,
const CellStatus,
- const boost::iterator_range<std::vector<char>::const_iterator>)>
+ const boost::iterator_range<std::vector<char>::const_iterator> &)>
&unpack_callback)
{
# ifdef DEBUG
template <int spacedim>
unsigned int
Triangulation<1, spacedim>::register_data_attach(
- const std::function<
- void(const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
- const typename dealii::Triangulation<1, spacedim>::CellStatus,
- std::vector<char> &)> & /*pack_callback*/)
+ const std::function<std::vector<char>(
+ const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
+ const typename dealii::Triangulation<1, spacedim>::CellStatus)>
+ & /*pack_callback*/)
{
Assert(false, ExcNotImplemented());
return 0;
const std::function<
void(const typename dealii::Triangulation<1, spacedim>::cell_iterator &,
const typename dealii::Triangulation<1, spacedim>::CellStatus,
- const boost::iterator_range<std::vector<char>::const_iterator>)>
+ const boost::iterator_range<std::vector<char>::const_iterator> &)>
& /*unpack_callback*/)
{
Assert(false, ExcNotImplemented());
if (global_max_particles_per_cell > 0)
{
- const std::function<
- void(const typename Triangulation<dim, spacedim>::cell_iterator &,
- const typename Triangulation<dim, spacedim>::CellStatus,
- std::vector<char> &)>
+ const std::function<std::vector<char>(
+ const typename Triangulation<dim, spacedim>::cell_iterator &,
+ const typename Triangulation<dim, spacedim>::CellStatus)>
callback_function =
std::bind(&ParticleHandler<dim, spacedim>::store_particles,
std::cref(*this),
std::placeholders::_1,
- std::placeholders::_2,
- std::placeholders::_3);
+ std::placeholders::_2);
handle =
non_const_triangulation->register_data_attach(callback_function);
// data correctly. Only do this if something was actually stored.
if (serialization && (global_max_particles_per_cell > 0))
{
- const std::function<
- void(const typename Triangulation<dim, spacedim>::cell_iterator &,
- const typename Triangulation<dim, spacedim>::CellStatus,
- std::vector<char> &)>
+ const std::function<std::vector<char>(
+ const typename Triangulation<dim, spacedim>::cell_iterator &,
+ const typename Triangulation<dim, spacedim>::CellStatus)>
callback_function =
std::bind(&ParticleHandler<dim, spacedim>::store_particles,
std::cref(*this),
std::placeholders::_1,
- std::placeholders::_2,
- std::placeholders::_3);
+ std::placeholders::_2);
handle =
non_const_triangulation->register_data_attach(callback_function);
// Check if something was stored and load it
if (handle != numbers::invalid_unsigned_int)
{
- const std::function<
- void(const typename Triangulation<dim, spacedim>::cell_iterator &,
- const typename Triangulation<dim, spacedim>::CellStatus,
- const boost::iterator_range<std::vector<char>::const_iterator>)>
+ const std::function<void(
+ const typename Triangulation<dim, spacedim>::cell_iterator &,
+ const typename Triangulation<dim, spacedim>::CellStatus,
+ const boost::iterator_range<std::vector<char>::const_iterator> &)>
callback_function =
std::bind(&ParticleHandler<dim, spacedim>::load_particles,
std::ref(*this),
template <int dim, int spacedim>
- void
+ std::vector<char>
ParticleHandler<dim, spacedim>::store_particles(
const typename Triangulation<dim, spacedim>::cell_iterator &cell,
- const typename Triangulation<dim, spacedim>::CellStatus status,
- std::vector<char> &data_vector) const
+ const typename Triangulation<dim, spacedim>::CellStatus status) const
{
// -------------------
// HOTFIX: resize memory and static_cast std::vector to void*
(size_per_particle * global_max_particles_per_cell) *
(serialization ? 1 : std::pow(2, dim));
- const size_t previous_size = data_vector.size();
- data_vector.resize(previous_size + transfer_size_per_cell);
- void *data = static_cast<void *>(data_vector.data() + previous_size);
+ std::vector<char> data_vector(transfer_size_per_cell);
+ void * data = static_cast<void *>(data_vector.data());
// --------------------
unsigned int n_particles(0);
}
else
Assert(false, ExcInternalError());
+
+ return data_vector;
}
template <int dim, int spacedim>
void
ParticleHandler<dim, spacedim>::load_particles(
- const typename Triangulation<dim, spacedim>::cell_iterator & cell,
- const typename Triangulation<dim, spacedim>::CellStatus status,
- const boost::iterator_range<std::vector<char>::const_iterator> data_range)
+ const typename Triangulation<dim, spacedim>::cell_iterator & cell,
+ const typename Triangulation<dim, spacedim>::CellStatus status,
+ const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
{
// -------------------
// HOTFIX: static_cast std::vector to void*
template <int dim>
-void
+std::vector<char>
pack_function(
const typename parallel::distributed::Triangulation<dim, dim>::cell_iterator
&cell,
const typename parallel::distributed::Triangulation<dim, dim>::CellStatus
- status,
- std::vector<char> &data)
+ status)
{
static int some_number = 0;
deallog << "packing cell " << cell->id() << " with data=" << some_number
Assert(!cell->has_children(), ExcInternalError());
}
- Utilities::pack(some_number, data, /*allow_compression=*/false);
-
- ++some_number;
+ return Utilities::pack(some_number++, /*allow_compression=*/false);
}
template <int dim>
const typename parallel::distributed::Triangulation<dim, dim>::cell_iterator
&cell,
const typename parallel::distributed::Triangulation<dim, dim>::CellStatus
- status,
- const boost::iterator_range<std::vector<char>::const_iterator> data_range)
+ status,
+ const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
{
const int intdata = Utilities::unpack<int>(data_range.begin(),
data_range.end(),
template <int dim>
-void
+std::vector<char>
pack_function(
const typename parallel::distributed::Triangulation<dim, dim>::cell_iterator
&cell,
const typename parallel::distributed::Triangulation<dim, dim>::CellStatus
- status,
- std::vector<char> &data)
+ status)
{
static int some_number = cell->index();
deallog << "packing cell " << cell->id() << " with data=" << some_number
Assert(!cell->has_children(), ExcInternalError());
}
- Utilities::pack(some_number, data, /*allow_compression=*/false);
-
- ++some_number;
+ return Utilities::pack(some_number++, /*allow_compression=*/false);
}
template <int dim>
const typename parallel::distributed::Triangulation<dim, dim>::cell_iterator
&cell,
const typename parallel::distributed::Triangulation<dim, dim>::CellStatus
- status,
- const boost::iterator_range<std::vector<char>::const_iterator> data_range)
+ status,
+ const boost::iterator_range<std::vector<char>::const_iterator> &data_range)
{
const int number = Utilities::unpack<int>(data_range.begin(),
data_range.end(),