#include <deal.II/particles/particle_accessor.h>
|
void * | write_particle_data_to_memory (void *data) const |
|
const void * | read_particle_data_from_memory (const void *data) |
|
void | set_location (const Point< spacedim > &new_location) |
|
const Point< spacedim > & | get_location () const |
|
void | set_reference_location (const Point< dim > &new_reference_location) |
|
const Point< dim > & | get_reference_location () const |
|
void | set_id (const types::particle_index &new_id) |
|
types::particle_index | get_id () const |
|
types::particle_index | get_local_index () const |
|
bool | has_properties () const |
|
void | set_properties (const std::vector< double > &new_properties) |
|
void | set_properties (const ArrayView< const double > &new_properties) |
|
void | set_properties (const Tensor< 1, dim > &new_properties) |
|
const ArrayView< double > | get_properties () |
|
const ArrayView< const double > | get_properties () const |
|
void | set_property_pool (PropertyPool< dim, spacedim > &property_pool) |
|
std::size_t | serialized_size_in_bytes () const |
|
const Triangulation< dim, spacedim >::cell_iterator & | get_surrounding_cell () const |
|
const Triangulation< dim, spacedim >::cell_iterator & | get_surrounding_cell (const Triangulation< dim, spacedim > &triangulation) const |
|
template<class Archive > |
void | save (Archive &ar, const unsigned int version) const |
|
template<class Archive > |
void | load (Archive &ar, const unsigned int version) |
|
template<class Archive > |
void | serialize (Archive &archive, const unsigned int version) |
|
void | next () |
|
void | prev () |
|
bool | operator!= (const ParticleAccessor< dim, spacedim > &other) const |
|
bool | operator== (const ParticleAccessor< dim, spacedim > &other) const |
|
IteratorState::IteratorStates | state () const |
|
template<
int dim,
int spacedim = dim>
class Particles::ParticleAccessor< dim, spacedim >
Accessor class used by ParticleIterator to access particle data.
Definition at line 45 of file particle_accessor.h.
◆ particle_container
template<
int dim,
int spacedim = dim>
◆ ParticleAccessor() [1/2]
template<
int dim,
int spacedim>
Construct an invalid accessor. Such an object is not usable.
Definition at line 540 of file particle_accessor.h.
◆ ParticleAccessor() [2/2]
template<
int dim,
int spacedim>
Particles::ParticleAccessor< dim, spacedim >::ParticleAccessor |
( |
const typename particle_container::iterator |
particles_in_cell, |
|
|
const PropertyPool< dim, spacedim > & |
property_pool, |
|
|
const unsigned int |
particle_index_within_cell |
|
) |
| |
|
inlineprivate |
Construct an accessor from a reference to a container, an iterator to the current cell, and the particle index within that cell. This constructor is private
so that it can only be accessed by friend classes.
Definition at line 549 of file particle_accessor.h.
◆ write_particle_data_to_memory()
template<
int dim,
int spacedim>
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 entry of the array to which the data should be written. This function is meant for serializing all particle properties and later de-serializing the properties by calling the appropriate constructor Particle(void *&data, PropertyPool *property_pool = nullptr);
- Parameters
-
[in] | data | The memory location to write particle data into. |
- Returns
- A pointer to the next byte after the array to which data has been written.
Definition at line 599 of file particle_accessor.h.
◆ read_particle_data_from_memory()
template<
int dim,
int spacedim>
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 which the data should be written. This function is meant for de-serializing the particle data without requiring that a new Particle class be built. This is used in the ParticleHandler to update the ghost particles without de-allocating and re-allocating memory.
- Parameters
-
[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. |
- Returns
- A pointer to the next byte after the array from which data has been read.
Definition at line 562 of file particle_accessor.h.
◆ set_location()
template<
int dim,
int spacedim>
Set the location of this particle. Note that this does not check whether this is a valid location in the simulation domain.
- Parameters
-
[in] | new_location | The new location for this particle. |
- Note
- In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as
const
objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.
Definition at line 633 of file particle_accessor.h.
◆ get_location()
template<
int dim,
int spacedim>
Get the location of this particle.
- Returns
- The location of this particle.
Definition at line 644 of file particle_accessor.h.
◆ set_reference_location()
template<
int dim,
int spacedim>
Set the reference location of this particle.
- Parameters
-
[in] | new_reference_location | The new reference location for this particle. |
- Note
- In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as
const
objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.
Definition at line 655 of file particle_accessor.h.
◆ get_reference_location()
template<
int dim,
int spacedim>
Return the reference location of this particle in its current cell.
Definition at line 667 of file particle_accessor.h.
◆ set_id()
template<
int dim,
int spacedim>
◆ get_id()
template<
int dim,
int spacedim>
◆ get_local_index()
template<
int dim,
int spacedim>
◆ has_properties()
template<
int dim,
int spacedim>
Return whether this particle has a valid property pool and a valid handle to properties.
Definition at line 711 of file particle_accessor.h.
◆ set_properties() [1/3]
template<
int dim,
int spacedim>
Set the properties of this particle.
- Parameters
-
[in] | new_properties | A vector containing the new properties for this particle. |
- Note
- In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as
const
objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.
Definition at line 729 of file particle_accessor.h.
◆ set_properties() [2/3]
template<
int dim,
int spacedim>
Set the properties of this particle.
- Parameters
-
[in] | new_properties | An ArrayView pointing to memory locations containing the new properties for this particle. |
- Note
- In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as
const
objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.
Definition at line 742 of file particle_accessor.h.
◆ set_properties() [3/3]
template<
int dim,
int spacedim>
Set the properties of this particle, assuming that the properties stored on this particle correspond to a rank-1 Tensor object. In particular, this means that the number of properties stored on the particle must equal dim
.
- Parameters
-
[in] | new_properties | A Tensor containing the new properties for this particle. |
- Note
- In parallel programs, the ParticleHandler class stores particles on both the locally owned cells, as well as on ghost cells. The particles on the latter are copies of particles owned on other processors, and should therefore be treated in the same way as ghost entries in vectors with ghost elements or ghost cells: In both cases, one should treat the ghost elements or cells as
const
objects that shouldn't be modified even if the objects allow for calls that modify properties. Rather, properties should only be modified on processors that actually own the particle.
Definition at line 769 of file particle_accessor.h.
◆ get_properties() [1/2]
template<
int dim,
int spacedim>
Get write-access to properties of this particle.
- Returns
- An ArrayView of the properties of this particle.
Definition at line 837 of file particle_accessor.h.
◆ get_properties() [2/2]
template<
int dim,
int spacedim>
Get read-access to properties of this particle.
- Returns
- An ArrayView of the properties of this particle.
Definition at line 788 of file particle_accessor.h.
◆ set_property_pool()
template<
int dim,
int spacedim>
Tell the particle where to store its properties (even if it does not own properties). Usually this is only done once per particle, but since the particle generator does not know about the properties we want to do it not at construction time. Another use for this function is after particle transfer to a new process.
- Deprecated:
- This function is only kept for backward compatibility and has no meaning any more. ParticleAccessors always use the property pool of the owning particle handler.
Definition at line 799 of file particle_accessor.h.
◆ serialized_size_in_bytes()
template<
int dim,
int spacedim>
Return the size in bytes this particle occupies if all of its data is serialized (i.e. the number of bytes that is written by the write_data function of this class).
Definition at line 848 of file particle_accessor.h.
◆ get_surrounding_cell() [1/2]
template<
int dim,
int spacedim>
Get a cell iterator to the cell surrounding the current particle. As particles are organized in the structure of a triangulation, but the triangulation itself is not stored in the particle this operation requires a reference to the triangulation.
Definition at line 810 of file particle_accessor.h.
◆ get_surrounding_cell() [2/2]
template<
int dim,
int spacedim>
◆ save()
template<
int dim,
int spacedim>
template<class Archive >
◆ load()
template<
int dim,
int spacedim>
template<class Archive >
Read the data of this object from a stream for the purpose of serialization using the BOOST serialization library. Note that in order to store the properties correctly, the property pool of this particle has to be known at the time of reading, i.e. set_property_pool() has to have been called, before this function is called.
Definition at line 484 of file particle_accessor.h.
◆ serialize()
template<
int dim,
int spacedim = dim>
template<class Archive >
◆ next()
template<
int dim,
int spacedim>
◆ prev()
template<
int dim,
int spacedim>
◆ operator!=()
template<
int dim,
int spacedim>
◆ operator==()
template<
int dim,
int spacedim>
◆ state()
template<
int dim,
int spacedim>
◆ get_handle() [1/2]
template<
int dim,
int spacedim>
Returns a reference to the current Particle. Because the internal structure may change this is not intended for public use and only a convenience function for internal purposes.
Definition at line 953 of file particle_accessor.h.
◆ get_handle() [2/2]
template<
int dim,
int spacedim>
Non-const
version of the function with the same name above.
Definition at line 944 of file particle_accessor.h.
◆ ParticleIterator
template<
int dim,
int spacedim = dim>
◆ ParticleHandler
template<
int dim,
int spacedim = dim>
◆ particles_in_cell
template<
int dim,
int spacedim = dim>
An iterator to the particles in the current cell within the particle_container object.
Definition at line 459 of file particle_accessor.h.
◆ property_pool
template<
int dim,
int spacedim = dim>
A pointer to the property pool that stores the actual particle data.
Definition at line 464 of file particle_accessor.h.
◆ particle_index_within_cell
template<
int dim,
int spacedim = dim>
The documentation for this class was generated from the following file: