* and later de-serializing the properties by calling the appropriate
* constructor Particle(void *&data, PropertyPool *property_pool = nullptr);
*
- * @param [in,out] data The memory location to write particle data
- * into. This pointer points to the begin of the memory, in which the
- * data will be written and it will be advanced by the serialized size
- * of this particle.
+ * @param [in] data The memory location to write particle data
+ * into.
+ *
+ * @return A pointer to the next byte after the array to which data has
+ * been written.
*/
- void
- write_data(void *&data) const;
+ void *
+ write_particle_data_to_memory(void *data) const;
/**
- * Update all of the data associated with a particle : id,
+ * Update all of the data associated with a particle: id,
* location, reference location and, if any, properties by using a
* data array. The array is expected to be large enough to take the data,
* and the void pointer should point to the first entry of the array to
* class be built. This is used in the ParticleHandler to update the
* ghost particles without de-allocating and re-allocating memory.
*
- * @param[in,out] data A pointer to a memory location from which
+ * @param[in] data A pointer to a memory location from which
* to read the information that completely describes a particle. This
- * class then de-serializes its data from this memory location and
- * advance the pointer accordingly.
+ * class then de-serializes its data from this memory location.
+ *
+ * @return A pointer to the next byte after the array from which data has
+ * been read.
*/
- void
- update_particle_data(const void *&data);
+ const void *
+ read_particle_data_from_memory(const void *data);
/**
* Set the location of this particle. Note that this does not check
{
public:
/**
- * Write particle data into a data array. The array is expected
- * to be large enough to take the data, and the void pointer should
- * point to the first element in which the data should be written. This
- * function is meant for serializing all particle properties and
- * afterwards de-serializing the properties by calling the appropriate
- * constructor Particle(void *&data, PropertyPool *property_pool = nullptr);
- *
- * @param [in,out] data The memory location to write particle data
- * into. This pointer points to the begin of the memory, in which the
- * data will be written and it will be advanced by the serialized size
- * of this particle.
+ * @copydoc Particle::write_particle_data_to_memory
*/
- void
- write_data(void *&data) const;
+ void *
+ write_particle_data_to_memory(void *data) const;
+
+
+ /**
+ * @copydoc Particle::read_particle_data_from_memory
+ */
+ const void *
+ read_particle_data_from_memory(const void *data);
/**
* Set the location of this particle. Note that this does not check
template <int dim, int spacedim>
- inline void
- ParticleAccessor<dim, spacedim>::write_data(void *&data) const
+ inline const void *
+ ParticleAccessor<dim, spacedim>::read_particle_data_from_memory(
+ const void *data)
+ {
+ Assert(particle != map->end(), ExcInternalError());
+
+ return particle->second.read_particle_data_from_memory(data);
+ }
+
+
+
+ template <int dim, int spacedim>
+ inline void *
+ ParticleAccessor<dim, spacedim>::write_particle_data_to_memory(
+ void *data) const
{
Assert(particle != map->end(), ExcInternalError());
- particle->second.write_data(data);
+ return particle->second.write_particle_data_to_memory(data);
}
template <int dim, int spacedim>
- void
- Particle<dim, spacedim>::write_data(void *&data) const
+ void *
+ Particle<dim, spacedim>::write_particle_data_to_memory(
+ void *data_pointer) const
{
- types::particle_index *id_data = static_cast<types::particle_index *>(data);
- *id_data = get_id();
+ types::particle_index *id_data =
+ static_cast<types::particle_index *>(data_pointer);
+ *id_data = get_id();
++id_data;
double *pdata = reinterpret_cast<double *>(id_data);
- // Write location data
+ // Write location
for (unsigned int i = 0; i < spacedim; ++i, ++pdata)
*pdata = get_location()[i];
- // Write reference location data
+ // Write reference location
for (unsigned int i = 0; i < dim; ++i, ++pdata)
*pdata = get_reference_location()[i];
- // Write property data
+ // Write properties
if (has_properties())
{
const ArrayView<double> particle_properties =
*pdata = particle_properties[i];
}
- data = static_cast<void *>(pdata);
+ return static_cast<void *>(pdata);
}
template <int dim, int spacedim>
- void
- Particle<dim, spacedim>::update_particle_data(const void *&data)
+ const void *
+ Particle<dim, spacedim>::read_particle_data_from_memory(
+ const void *data_pointer)
{
const types::particle_index *id_data =
- static_cast<const types::particle_index *>(data);
+ static_cast<const types::particle_index *>(data_pointer);
set_id(*id_data++);
const double *pdata = reinterpret_cast<const double *>(id_data);
particle_properties[i] = *pdata++;
}
- data = static_cast<const void *>(pdata);
+ return static_cast<const void *>(pdata);
}
for (const auto &particle : particles)
{
- particle.write_data(current_data);
+ current_data = particle.write_particle_data_to_memory(current_data);
}
return buffer;
}
+
+
template <int dim, int spacedim>
std::vector<Particle<dim, spacedim>>
unpack_particles(
memcpy(data, &cellid, cellid_size);
data = static_cast<char *>(data) + cellid_size;
- particles_to_send.at(neighbors[i])[j]->write_data(data);
+ data = particles_to_send.at(neighbors[i])[j]
+ ->write_particle_data_to_memory(data);
if (store_callback)
data =
store_callback(particles_to_send.at(neighbors[i])[j], data);
for (const auto i : neighbors)
for (const auto &p : particles_to_send.at(i))
{
- p->write_data(data);
+ data = p->write_particle_data_to_memory(data);
if (store_callback)
data = store_callback(p, data);
}
{
// Update particle data using previously allocated memory space
// for efficiency reasons
- recv_particle->second.update_particle_data(recv_data_it);
+ recv_data_it =
+ recv_particle->second.read_particle_data_from_memory(recv_data_it);
if (load_callback)
recv_data_it =
return pack_particles(stored_particles_on_cell);
}
+
+
template <int dim, int spacedim>
void
ParticleHandler<dim, spacedim>::load_particles(
std::vector<char> data(particle.serialized_size_in_bytes());
void * write_pointer = static_cast<void *>(&data.front());
- particle.write_data(write_pointer);
+ write_pointer = particle.write_particle_data_to_memory(write_pointer);
const void *read_pointer = static_cast<const void *>(&data.front());
const Particles::Particle<dim, spacedim> new_particle(read_pointer);
std::vector<char> data(particle.serialized_size_in_bytes());
void * write_pointer = static_cast<void *>(&data.front());
- particle.write_data(write_pointer);
+ write_pointer = particle.write_particle_data_to_memory(write_pointer);
const void *read_pointer = static_cast<const void *>(&data.front());
const Particles::Particle<dim, spacedim> new_particle(read_pointer, &pool);