Reference documentation for deal.II version GIT relicensing-249-g48dc7357c7 2024-03-29 12:30:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Classes | Public Types | Public Member Functions | Private Member Functions | Private Attributes | Friends | List of all members
Particles::ParticleAccessor< dim, spacedim > Class Template Reference

#include <deal.II/particles/particle_accessor.h>

Inheritance diagram for Particles::ParticleAccessor< dim, spacedim >:
Inheritance graph
[legend]

Classes

struct  ParticlesInCell
 

Public Types

using particle_container = std::list< ParticlesInCell >
 

Public Member Functions

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
 
Point< spacedim > & get_location ()
 
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)
 
ArrayView< double > get_properties ()
 
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
 

Private Member Functions

 ParticleAccessor ()
 
 ParticleAccessor (const typename particle_container::iterator particles_in_cell, const PropertyPool< dim, spacedim > &property_pool, const unsigned int particle_index_within_cell)
 
const PropertyPool< dim, spacedim >::Handle & get_handle () const
 
PropertyPool< dim, spacedim >::Handle & get_handle ()
 

Private Attributes

particle_container::iterator particles_in_cell
 
PropertyPool< dim, spacedim > * property_pool
 
unsigned int particle_index_within_cell
 

Friends

template<int , int >
class ParticleIterator
 
template<int , int >
class ParticleHandler
 

Detailed Description

template<int dim, int spacedim = dim>
class Particles::ParticleAccessor< dim, spacedim >

Accessor class used by ParticleIterator to access particle data.

Definition at line 44 of file particle_accessor.h.

Member Typedef Documentation

◆ particle_container

template<int dim, int spacedim = dim>
using Particles::ParticleAccessor< dim, spacedim >::particle_container = std::list<ParticlesInCell>

A type for the storage container for particles.

Definition at line 112 of file particle_accessor.h.

Constructor & Destructor Documentation

◆ ParticleAccessor() [1/2]

template<int dim, int spacedim>
Particles::ParticleAccessor< dim, spacedim >::ParticleAccessor ( )
inlineprivate

Construct an invalid accessor. Such an object is not usable.

Definition at line 563 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 572 of file particle_accessor.h.

Member Function Documentation

◆ write_particle_data_to_memory()

template<int dim, int spacedim>
void * Particles::ParticleAccessor< dim, spacedim >::write_particle_data_to_memory ( void *  data) const
inline

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]dataThe 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 622 of file particle_accessor.h.

◆ read_particle_data_from_memory()

template<int dim, int spacedim>
const void * Particles::ParticleAccessor< dim, spacedim >::read_particle_data_from_memory ( const void *  data)
inline

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]dataA 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 585 of file particle_accessor.h.

◆ set_location()

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::set_location ( const Point< spacedim > &  new_location)
inline

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_locationThe 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 656 of file particle_accessor.h.

◆ get_location() [1/2]

template<int dim, int spacedim>
const Point< spacedim > & Particles::ParticleAccessor< dim, spacedim >::get_location ( ) const
inline

Get read-access to the location of this particle.

Returns
The location of this particle.

Definition at line 667 of file particle_accessor.h.

◆ get_location() [2/2]

template<int dim, int spacedim>
Point< spacedim > & Particles::ParticleAccessor< dim, spacedim >::get_location ( )
inline

Get read- and write-access to the location of this particle. Note that changing the location does not check whether this is a valid location in the simulation domain.

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.
Returns
The location of this particle.

Definition at line 678 of file particle_accessor.h.

◆ set_reference_location()

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::set_reference_location ( const Point< dim > &  new_reference_location)
inline

Set the reference location of this particle.

Parameters
[in]new_reference_locationThe 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 689 of file particle_accessor.h.

◆ get_reference_location()

template<int dim, int spacedim>
const Point< dim > & Particles::ParticleAccessor< dim, spacedim >::get_reference_location ( ) const
inline

Return the reference location of this particle in its current cell.

Definition at line 701 of file particle_accessor.h.

◆ set_id()

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::set_id ( const types::particle_index new_id)
inline

Set the ID number of this particle.

Definition at line 734 of file particle_accessor.h.

◆ get_id()

template<int dim, int spacedim>
types::particle_index Particles::ParticleAccessor< dim, spacedim >::get_id ( ) const
inline

Return the ID number of this particle.

Definition at line 712 of file particle_accessor.h.

◆ get_local_index()

template<int dim, int spacedim>
types::particle_index Particles::ParticleAccessor< dim, spacedim >::get_local_index ( ) const
inline

Return a particle ID number local to each MPI process. This number enables the direct array access (similar to LinearAlgebra::distributed::Vector::local_element()) to quantities used for local computations. Use ParticleHandler::get_max_local_particle_index() to query suitable array sizes.

Note
The number returned by this function is not stable between calls of ParticleHandler::sort_particles_into_subdomains_and_cells() and therefore it cannot be used to track quantities across those calls. Furthermore, the numbers returned by this function are typically not refreshed in case individual particles are removed from the underlying ParticleHandler object, which means that the returned indices of all particles on an MPI process typically do not form a contiguous interval of numbers. Use particle properties for storing persistent information and checkpointing the actual data.

Definition at line 723 of file particle_accessor.h.

◆ has_properties()

template<int dim, int spacedim>
bool Particles::ParticleAccessor< dim, spacedim >::has_properties ( ) const
inline

Return whether this particle has a valid property pool and a valid handle to properties.

Definition at line 745 of file particle_accessor.h.

◆ set_properties() [1/3]

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::set_properties ( const std::vector< double > &  new_properties)
inline

Set the properties of this particle.

Parameters
[in]new_propertiesA 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 763 of file particle_accessor.h.

◆ set_properties() [2/3]

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::set_properties ( const ArrayView< const double > &  new_properties)
inline

Set the properties of this particle.

Parameters
[in]new_propertiesAn 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 776 of file particle_accessor.h.

◆ set_properties() [3/3]

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::set_properties ( const Tensor< 1, dim > &  new_properties)
inline

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_propertiesA 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 803 of file particle_accessor.h.

◆ get_properties() [1/2]

template<int dim, int spacedim>
ArrayView< double > Particles::ParticleAccessor< dim, spacedim >::get_properties ( )
inline

Get write-access to properties of this particle.

Returns
An ArrayView of the properties of this particle.

Definition at line 871 of file particle_accessor.h.

◆ get_properties() [2/2]

template<int dim, int spacedim>
ArrayView< const double > Particles::ParticleAccessor< dim, spacedim >::get_properties ( ) const
inline

Get read-access to properties of this particle.

Returns
An ArrayView of the properties of this particle.

Definition at line 822 of file particle_accessor.h.

◆ set_property_pool()

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::set_property_pool ( PropertyPool< dim, spacedim > &  property_pool)
inline

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 833 of file particle_accessor.h.

◆ serialized_size_in_bytes()

template<int dim, int spacedim>
std::size_t Particles::ParticleAccessor< dim, spacedim >::serialized_size_in_bytes ( ) const
inline

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 882 of file particle_accessor.h.

◆ get_surrounding_cell() [1/2]

template<int dim, int spacedim>
const Triangulation< dim, spacedim >::cell_iterator & Particles::ParticleAccessor< dim, spacedim >::get_surrounding_cell ( ) const
inline

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 844 of file particle_accessor.h.

◆ get_surrounding_cell() [2/2]

template<int dim, int spacedim>
const Triangulation< dim, spacedim >::cell_iterator & Particles::ParticleAccessor< dim, spacedim >::get_surrounding_cell ( const Triangulation< dim, spacedim > &  triangulation) const
inline
Deprecated:
Deprecated version of the function with the same name above.

Definition at line 857 of file particle_accessor.h.

◆ save()

template<int dim, int spacedim>
template<class Archive >
void Particles::ParticleAccessor< dim, spacedim >::save ( Archive &  ar,
const unsigned int  version 
) const
inline

Write the data of this object to a stream for the purpose of serialization using the BOOST serialization library.

Definition at line 541 of file particle_accessor.h.

◆ load()

template<int dim, int spacedim>
template<class Archive >
void Particles::ParticleAccessor< dim, spacedim >::load ( Archive &  ar,
const unsigned int  version 
)
inline

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 507 of file particle_accessor.h.

◆ serialize()

template<int dim, int spacedim = dim>
template<class Archive >
void Particles::ParticleAccessor< dim, spacedim >::serialize ( Archive &  archive,
const unsigned int  version 
)

Write and read the data of this object from a stream for the purpose of serialization using the BOOST serialization library.

◆ next()

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::next ( )
inline

Advance the ParticleAccessor to the next particle.

Definition at line 900 of file particle_accessor.h.

◆ prev()

template<int dim, int spacedim>
void Particles::ParticleAccessor< dim, spacedim >::prev ( )
inline

Move the ParticleAccessor to the previous particle.

Definition at line 917 of file particle_accessor.h.

◆ operator!=()

template<int dim, int spacedim>
bool Particles::ParticleAccessor< dim, spacedim >::operator!= ( const ParticleAccessor< dim, spacedim > &  other) const
inline

Inequality operator.

Definition at line 936 of file particle_accessor.h.

◆ operator==()

template<int dim, int spacedim>
bool Particles::ParticleAccessor< dim, spacedim >::operator== ( const ParticleAccessor< dim, spacedim > &  other) const
inline

Equality operator.

Definition at line 946 of file particle_accessor.h.

◆ state()

template<int dim, int spacedim>
IteratorState::IteratorStates Particles::ParticleAccessor< dim, spacedim >::state ( ) const
inline

Return the state of the accessor.

Definition at line 958 of file particle_accessor.h.

◆ get_handle() [1/2]

template<int dim, int spacedim>
const PropertyPool< dim, spacedim >::Handle & Particles::ParticleAccessor< dim, spacedim >::get_handle ( ) const
inlineprivate

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 987 of file particle_accessor.h.

◆ get_handle() [2/2]

template<int dim, int spacedim>
PropertyPool< dim, spacedim >::Handle & Particles::ParticleAccessor< dim, spacedim >::get_handle ( )
inlineprivate

Non-const version of the function with the same name above.

Definition at line 978 of file particle_accessor.h.

Friends And Related Symbol Documentation

◆ ParticleIterator

template<int dim, int spacedim = dim>
template<int , int >
friend class ParticleIterator
friend

Definition at line 497 of file particle_accessor.h.

◆ ParticleHandler

template<int dim, int spacedim = dim>
template<int , int >
friend class ParticleHandler
friend

Definition at line 499 of file particle_accessor.h.

Member Data Documentation

◆ particles_in_cell

template<int dim, int spacedim = dim>
particle_container::iterator Particles::ParticleAccessor< dim, spacedim >::particles_in_cell
private

An iterator to the particles in the current cell within the particle_container object.

Definition at line 482 of file particle_accessor.h.

◆ property_pool

template<int dim, int spacedim = dim>
PropertyPool<dim, spacedim>* Particles::ParticleAccessor< dim, spacedim >::property_pool
private

A pointer to the property pool that stores the actual particle data.

Definition at line 487 of file particle_accessor.h.

◆ particle_index_within_cell

template<int dim, int spacedim = dim>
unsigned int Particles::ParticleAccessor< dim, spacedim >::particle_index_within_cell
private

Local index of the particle within its current cell.

Definition at line 492 of file particle_accessor.h.


The documentation for this class was generated from the following file: