FullMatrix<double> matrix_quadrature;
/**
- * The offset that the parallel::distributed::Triangulation has assigned
- * to this object starting at which we are allowed to read.
+ * The handle that the parallel::distributed::Triangulation has assigned
+ * to this object while registering the pack_callback function.
*/
- unsigned int offset;
+ unsigned int handle;
/**
* A pointer to the CellDataStorage class whose data will be transferred.
n_q_points(rhs_quadrature.size()),
project_to_fe_matrix(projection_fe->dofs_per_cell,n_q_points),
project_to_qp_matrix(n_q_points,projection_fe->dofs_per_cell),
- offset(0),
+ handle(numbers::invalid_unsigned_int),
data_storage(nullptr),
triangulation(nullptr)
{
matrix_quadrature.reinit(n_q_points,number_of_values);
data_size_in_bytes = sizeof(double) * dofs_per_cell * number_of_values;
- offset = triangulation->register_data_attach(data_size_in_bytes,
+ handle = triangulation->register_data_attach(data_size_in_bytes,
std::bind(&ContinuousQuadratureDataTransfer<dim,DataType>::pack_function,
this,
std::placeholders::_1,
template <int dim, typename DataType>
void ContinuousQuadratureDataTransfer<dim,DataType>::interpolate ()
{
- triangulation->notify_ready_to_unpack(offset,
+ triangulation->notify_ready_to_unpack(handle,
std::bind(&ContinuousQuadratureDataTransfer<dim,DataType>::unpack_function,
this,
std::placeholders::_1,
std::vector<const VectorType *> input_vectors;
/**
- * The offset that the Triangulation has assigned to this object
- * starting at which we are allowed to write.
+ * The handle that the Triangulation has assigned to this object
+ * with which we can access our memory offset and our pack function.
*/
- unsigned int offset;
+ unsigned int handle;
/**
* A callback function used to pack the data on the current mesh into
/**
* List of callback functions registered by register_data_attach() that
- * are going to be called for packing data. The key of this map
- * variable is the offset at which each callback is allowed to
- * write into the per-cell buffer (counted in bytes) whereas the
- * value of the map is the callback function object.
+ * are going to be called for packing data. The callback functions objects
+ * will be stored together with their offset at which each callback is
+ * allowed to write into the per-cell buffer (counted in bytes). These
+ * pairs will be stored in the order on how they have been registered.
*/
- std::map<unsigned int, pack_callback_t> pack_callbacks;
+ std::vector< std::pair<unsigned int, pack_callback_t> > pack_callbacks;
};
CellAttachedData cell_attached_data;
/**
* This variable is set by the register_store_callback_function()
* function and used by the register_load_callback_function() function
- * to check where the particle data was stored.
+ * to check where the particle data was registered in the corresponding
+ * triangulation object.
*/
- unsigned int data_offset;
+ unsigned int handle;
#ifdef DEAL_II_WITH_MPI
/**
SolutionTransfer<dim, VectorType, DoFHandlerType>::SolutionTransfer (const DoFHandlerType &dof)
:
dof_handler(&dof, typeid(*this).name()),
- offset (numbers::invalid_unsigned_int)
+ handle (numbers::invalid_unsigned_int)
{
Assert ((dynamic_cast<const parallel::distributed::Triangulation<dim,DoFHandlerType::space_dimension>*>
(&dof_handler->get_triangulation()) != nullptr),
(&dof_handler->get_triangulation())));
Assert (tria != nullptr, ExcInternalError());
- offset
+ handle
= tria->register_data_attach(size,
std::bind(&SolutionTransfer<dim, VectorType,
DoFHandlerType>::pack_callback,
(&dof_handler->get_triangulation())));
Assert (tria != nullptr, ExcInternalError());
- tria->notify_ready_to_unpack(offset,
+ tria->notify_ready_to_unpack(handle,
std::bind(&SolutionTransfer<dim, VectorType,
DoFHandlerType>::unpack_callback,
this,
attach_mesh_data_recursively (const typename internal::p4est::types<dim>::tree &tree,
const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
const typename internal::p4est::types<dim>::quadrant &p4est_cell,
- const std::map<unsigned int, typename std::function<
+ const std::vector< std::pair<unsigned int, typename std::function<
void(typename parallel::distributed::Triangulation<dim,spacedim>::cell_iterator,
typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus,
void *)
- > > &pack_callbacks)
+ > > > &pack_callbacks)
{
int idx = sc_array_bsearch(const_cast<sc_array_t *>(&tree.quadrants),
&p4est_cell,
const typename Triangulation<dim,spacedim>::cell_iterator &dealii_cell,
const typename Triangulation<dim,spacedim>::cell_iterator &parent_cell,
const typename internal::p4est::types<dim>::quadrant &p4est_cell,
- const unsigned int handle,
+ const unsigned int offset,
const typename std::function<
void(typename parallel::distributed::Triangulation<dim,spacedim>::cell_iterator, typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus, void *)
> &unpack_callback)
dealii_cell->child(c),
dealii_cell,
p4est_child[c],
- handle,
+ offset,
unpack_callback);
}
}
sc_array_index (const_cast<sc_array_t *>(&tree.quadrants), idx)
);
- void *ptr = static_cast<char *>(q->p.user_data) + handle;
+ void *ptr = static_cast<char *>(q->p.user_data) + offset;
typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus
status = * static_cast<
typename parallel::distributed::Triangulation<dim,spacedim>::CellStatus *
const CellStatus,
void *)> &pack_callback)
{
+#ifdef DEBUG
Assert(size>0, ExcMessage("register_data_attach(), size==0"));
- Assert(cell_attached_data.pack_callbacks.size()==cell_attached_data.n_attached_datas,
+
+ unsigned int remaining_callback_counter = 0;
+ for (auto it : cell_attached_data.pack_callbacks)
+ if (std::get<0>(it) != numbers::invalid_unsigned_int)
+ ++remaining_callback_counter;
+
+ Assert(remaining_callback_counter==cell_attached_data.n_attached_datas,
ExcMessage("register_data_attach(), not all data has been unpacked last time?"));
+#endif
+
+ // add new pair with offset and pack_callback
+ cell_attached_data.pack_callbacks.push_back(
+ std::make_pair(cell_attached_data.cumulative_fixed_data_size + sizeof(CellStatus),
+ pack_callback)
+ );
- const unsigned int current_buffer_positition = cell_attached_data.cumulative_fixed_data_size+sizeof(CellStatus);
+ // increase counters
++cell_attached_data.n_attached_datas;
cell_attached_data.cumulative_fixed_data_size += size;
- cell_attached_data.pack_callbacks[current_buffer_positition] = pack_callback;
- // use the offset as the handle that callers receive and later hand back
- // to notify_ready_to_unpack()
- return current_buffer_positition;
+ // return the position of the callback in the register vector as a handle
+ return cell_attached_data.pack_callbacks.size() - 1;
}
const CellStatus,
const void *)> &unpack_callback)
{
- Assert (handle >= sizeof(CellStatus),
- ExcMessage ("Invalid offset."));
- Assert (handle < cell_attached_data.cumulative_fixed_data_size + sizeof(CellStatus),
- ExcMessage ("Invalid offset."));
+ Assert (handle < cell_attached_data.pack_callbacks.size(),
+ ExcMessage ("Invalid handle."));
Assert (cell_attached_data.n_attached_datas > 0,
ExcMessage ("The notify_ready_to_unpack() has already been called "
"once for each registered callback."));
+ const unsigned int offset = std::get<0> (cell_attached_data.pack_callbacks[handle]);
+ // check if unpack function has been previously called
+ Assert (offset!=dealii::numbers::invalid_unsigned_int, ExcInternalError());
+
// Recurse over p4est and hand the caller the data back
for (typename Triangulation<dim, spacedim>::cell_iterator
cell = this->begin (0);
cell,
cell,
p4est_coarse_cell,
- handle,
+ offset,
unpack_callback);
}
--cell_attached_data.n_attached_datas;
if (cell_attached_data.n_attached_deserialize > 0)
- {
- --cell_attached_data.n_attached_deserialize;
-
- // Remove the entry that corresponds to the handle that was
- // passed to the current function. Double check that such an
- // entry really existed.
- const unsigned int n_elements_removed
- = cell_attached_data.pack_callbacks.erase (handle);
- (void)n_elements_removed;
- Assert (n_elements_removed == 1,
- ExcInternalError());
- }
+ --cell_attached_data.n_attached_deserialize;
+
+ // deregister corresponding pack_callback function
+ // first reset with invalid entries to preserve ambiguity of handles,
+ // then free memory when all were unpacked (see below)
+ std::get<0> (cell_attached_data.pack_callbacks[handle]) // offset value
+ = dealii::numbers::invalid_unsigned_int;
+ std::get<1> (cell_attached_data.pack_callbacks[handle]) // pack_callback function
+ = typename CellAttachedData::pack_callback_t ();
// important: only remove data if we are not in the deserialization
// process. There, each SolutionTransfer registers and unpacks before
size_callback(),
store_callback(),
load_callback(),
- data_offset(numbers::invalid_unsigned_int)
+ handle(numbers::invalid_unsigned_int)
{}
size_callback(),
store_callback(),
load_callback(),
- data_offset(numbers::invalid_unsigned_int)
+ handle(numbers::invalid_unsigned_int)
{}
:
std::pow(2,dim));
- data_offset = non_const_triangulation->register_data_attach(transfer_size_per_cell,callback_function);
+ handle = non_const_triangulation->register_data_attach(transfer_size_per_cell,callback_function);
}
}
// space for the data in case a cell is coarsened
const std::size_t transfer_size_per_cell = sizeof (unsigned int) +
(size_per_particle * global_max_particles_per_cell);
- data_offset = non_const_triangulation->register_data_attach(transfer_size_per_cell,callback_function);
+ handle = non_const_triangulation->register_data_attach(transfer_size_per_cell,callback_function);
}
// Check if something was stored and load it
- if (data_offset != numbers::invalid_unsigned_int)
+ if (handle != numbers::invalid_unsigned_int)
{
const std::function<void(const typename Triangulation<dim,spacedim>::cell_iterator &,
const typename Triangulation<dim,spacedim>::CellStatus,
std::placeholders::_2,
std::placeholders::_3);
- non_const_triangulation->notify_ready_to_unpack(data_offset,callback_function);
+ non_const_triangulation->notify_ready_to_unpack(handle,callback_function);
- // Reset offset and update global number of particles. The number
+ // Reset handle and update global number of particles. The number
// can change because of discarded or newly generated particles
- data_offset = numbers::invalid_unsigned_int;
+ handle = numbers::invalid_unsigned_int;
update_cached_numbers();
}
}
}
}
- unsigned int offset = tr.register_data_attach(sizeof(int), pack_function<dim>);
+ unsigned int handle = tr.register_data_attach(sizeof(int), pack_function<dim>);
- deallog << "offset=" << offset << std::endl;
+ deallog << "handle=" << handle << std::endl;
tr.execute_coarsening_and_refinement();
deallog << "cells after: " << tr.n_global_active_cells() << std::endl;
deallog << "myid=" << myid << " cellid=" << cell->id() << std::endl;
}*/
- tr.notify_ready_to_unpack(offset, unpack_function<dim>);
+ tr.notify_ready_to_unpack(handle, unpack_function<dim>);
const unsigned int checksum = tr.get_checksum ();
deallog << "Checksum: "
DEAL:0::myid=0 cellid=3_1:1 coarsening
DEAL:0::myid=0 cellid=3_1:2 coarsening
DEAL:0::myid=0 cellid=3_1:3 coarsening
-DEAL:0::offset=4
+DEAL:0::handle=0
DEAL:0::packing cell 0_1:0 with data=0 status=REFINE
DEAL:0::packing cell 0_1:1 with data=1 status=PERSIST
DEAL:0::packing cell 0_1:2 with data=2 status=PERSIST
DEAL:0::myid=0 cellid=0_1:1
DEAL:0::myid=0 cellid=0_1:2
DEAL:0::myid=0 cellid=0_1:3
-DEAL:0::offset=4
+DEAL:0::handle=0
DEAL:0::packing cell 0_1:0 with data=0 status=REFINE
DEAL:0::packing cell 0_1:1 with data=1 status=PERSIST
DEAL:0::packing cell 0_1:2 with data=2 status=PERSIST
DEAL:1::myid=1 cellid=2_1:1
DEAL:1::myid=1 cellid=2_1:2
DEAL:1::myid=1 cellid=2_1:3
-DEAL:1::offset=4
+DEAL:1::handle=0
DEAL:1::packing cell 1_1:0 with data=0 status=PERSIST
DEAL:1::packing cell 1_1:1 with data=1 status=PERSIST
DEAL:1::packing cell 1_1:2 with data=2 status=PERSIST
DEAL:2::myid=2 cellid=3_1:1 coarsening
DEAL:2::myid=2 cellid=3_1:2 coarsening
DEAL:2::myid=2 cellid=3_1:3 coarsening
-DEAL:2::offset=4
+DEAL:2::handle=0
DEAL:2::packing cell 3_0: with data=0 status=COARSEN
DEAL:2::cells after: 16
DEAL:2::unpacking cell 2_1:0 with data=4 status=PERSIST