Reference documentation for deal.II version 9.0.0
|
#include <deal.II/particles/particle_handler.h>
Public Types | |
typedef ParticleIterator< dim, spacedim > | particle_iterator |
typedef boost::iterator_range< particle_iterator > | particle_iterator_range |
Public Member Functions | |
ParticleHandler () | |
ParticleHandler (const parallel::distributed::Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const unsigned int n_properties=0) | |
~ParticleHandler () | |
void | initialize (const parallel::distributed::Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const unsigned int n_properties=0) |
void | clear () |
void | clear_particles () |
void | update_cached_numbers () |
particle_iterator | begin () const |
particle_iterator | begin () |
particle_iterator | end () const |
particle_iterator | end () |
particle_iterator | begin_ghost () const |
particle_iterator | begin_ghost () |
particle_iterator | end_ghost () const |
particle_iterator | end_ghost () |
particle_iterator_range | particles_in_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) |
particle_iterator_range | particles_in_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const |
void | remove_particle (const particle_iterator &particle) |
particle_iterator | insert_particle (const Particle< dim, spacedim > &particle, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) |
void | insert_particles (const std::multimap< typename Triangulation< dim, spacedim >::active_cell_iterator, Particle< dim, spacedim > > &particles) |
void | insert_particles (const std::vector< Point< spacedim > > &positions) |
void | register_additional_store_load_functions (const std::function< std::size_t()> &size_callback, const std::function< void *(const particle_iterator &, void *)> &store_callback, const std::function< const void *(const particle_iterator &, const void *)> &load_callback) |
types::particle_index | n_global_particles () const |
types::particle_index | n_global_max_particles_per_cell () const |
types::particle_index | n_locally_owned_particles () const |
types::particle_index | get_next_free_particle_index () const |
unsigned int | n_properties_per_particle () const |
PropertyPool & | get_property_pool () const |
unsigned int | n_particles_in_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const |
void | sort_particles_into_subdomains_and_cells () |
void | exchange_ghost_particles () |
void | register_store_callback_function (const bool serialization) |
void | register_load_callback_function (const bool serialization) |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Public Member Functions inherited from Subscriptor | |
Subscriptor () | |
Subscriptor (const Subscriptor &) | |
Subscriptor (Subscriptor &&) noexcept | |
virtual | ~Subscriptor () |
Subscriptor & | operator= (const Subscriptor &) |
Subscriptor & | operator= (Subscriptor &&) noexcept |
void | subscribe (const char *identifier=nullptr) const |
void | unsubscribe (const char *identifier=nullptr) const |
unsigned int | n_subscriptions () const |
void | list_subscribers () const |
template<class Archive > | |
void | serialize (Archive &ar, const unsigned int version) |
Private Member Functions | |
void | send_recv_particles (const std::map< types::subdomain_id, std::vector< particle_iterator > > &particles_to_send, std::multimap< internal::LevelInd, Particle< dim, spacedim > > &received_particles, const std::map< types::subdomain_id, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > > &new_cells_for_particles=std::map< types::subdomain_id, std::vector< typename Triangulation< dim, spacedim >::active_cell_iterator > >()) |
void | store_particles (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status, void *data) const |
void | load_particles (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status, const void *data) |
Private Attributes | |
SmartPointer< const parallel::distributed::Triangulation< dim, spacedim >, ParticleHandler< dim, spacedim > > | triangulation |
SmartPointer< const Mapping< dim, spacedim >, ParticleHandler< dim, spacedim > > | mapping |
std::multimap< internal::LevelInd, Particle< dim, spacedim > > | particles |
std::multimap< internal::LevelInd, Particle< dim, spacedim > > | ghost_particles |
types::particle_index | global_number_of_particles |
unsigned int | global_max_particles_per_cell |
types::particle_index | next_free_particle_index |
std::unique_ptr< PropertyPool > | property_pool |
std::function< std::size_t()> | size_callback |
std::function< void *(const particle_iterator &, void *)> | store_callback |
std::function< const void *(const particle_iterator &, const void *)> | load_callback |
unsigned int | data_offset |
Additional Inherited Members | |
Static Public Member Functions inherited from Subscriptor | |
static ::ExceptionBase & | ExcInUse (int arg1, std::string arg2, std::string arg3) |
static ::ExceptionBase & | ExcNoSubscriber (std::string arg1, std::string arg2) |
This class manages the storage and handling of particles. It provides the data structures necessary to store particles efficiently, accessor functions to iterate over particles and find particles, and algorithms to distribute particles in parallel domains. Note that the class is designed in a similar way as the triangulation class. In particular, we call particles in the domain of the local process local particles, and particles that belong to neighbor processes and live in the ghost cells around the locally owned domain "ghost particles".
Definition at line 28 of file particle_accessor.h.
typedef ParticleIterator<dim,spacedim> Particles::ParticleHandler< dim, spacedim >::particle_iterator |
A type that can be used to iterate over all particles in the domain.
Definition at line 58 of file particle_handler.h.
typedef boost::iterator_range<particle_iterator> Particles::ParticleHandler< dim, spacedim >::particle_iterator_range |
A type that represents a range of particles.
Definition at line 63 of file particle_handler.h.
Particles::ParticleHandler< dim, spacedim >::ParticleHandler | ( | ) |
Default constructor.
Definition at line 31 of file particle_handler.cc.
Particles::ParticleHandler< dim, spacedim >::ParticleHandler | ( | const parallel::distributed::Triangulation< dim, spacedim > & | tria, |
const Mapping< dim, spacedim > & | mapping, | ||
const unsigned int | n_properties = 0 |
||
) |
Constructor that initializes the particle handler with a given triangulation and mapping. Since particles are stored in respect to their surrounding cells this information is necessary to correctly organize the particle collection. This constructor is equivalent to calling the default constructor and the initialize function.
Definition at line 49 of file particle_handler.cc.
Particles::ParticleHandler< dim, spacedim >::~ParticleHandler | ( | ) |
Destructor.
Definition at line 70 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::initialize | ( | const parallel::distributed::Triangulation< dim, spacedim > & | tria, |
const Mapping< dim, spacedim > & | mapping, | ||
const unsigned int | n_properties = 0 |
||
) |
Initialize the particle handler. This function does not clear the internal data structures, it just sets the connections to the MPI communicator and the triangulation.
Definition at line 77 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::clear | ( | ) |
Clear all particle related data.
Definition at line 92 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::clear_particles | ( | ) |
Only clear particle data, but keep cache information about number of particles. This is useful during reorganization of particle data between processes.
Definition at line 104 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::update_cached_numbers | ( | ) |
Update all internally cached numbers. Note that all functions that modify internal data structures and act on multiple particles will call this function automatically (e.g. insert_particles), while functions that act on single particles will not call this function (e.g. insert_particle). This is done because the update is expensive compared to single operations.
Definition at line 113 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::begin | ( | ) | const |
Return an iterator to the first particle.
Definition at line 146 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::begin | ( | ) |
Return an iterator to the first particle.
Definition at line 155 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::end | ( | ) | const |
Return an iterator past the end of the particles.
Definition at line 164 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::end | ( | ) |
Return an iterator past the end of the particles.
Definition at line 173 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::begin_ghost | ( | ) | const |
Return an iterator to the first ghost particle.
Definition at line 182 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::begin_ghost | ( | ) |
Return an iterator to the first ghost particle.
Definition at line 191 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::end_ghost | ( | ) | const |
Return an iterator past the end of the ghost particles.
Definition at line 200 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::end_ghost | ( | ) |
Return an iterator past the end of the ghost particles.
Definition at line 209 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator_range Particles::ParticleHandler< dim, spacedim >::particles_in_cell | ( | const typename Triangulation< dim, spacedim >::active_cell_iterator & | cell | ) |
Return a pair of particle iterators that mark the begin and end of the particles in a particular cell. The last iterator is the first particle that is no longer in the cell.
Definition at line 227 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator_range Particles::ParticleHandler< dim, spacedim >::particles_in_cell | ( | const typename Triangulation< dim, spacedim >::active_cell_iterator & | cell | ) | const |
Return a pair of particle iterators that mark the begin and end of the particles in a particular cell. The last iterator is the first particle that is no longer in the cell.
Definition at line 218 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::remove_particle | ( | const particle_iterator & | particle | ) |
Remove a particle pointed to by the iterator.
Definition at line 247 of file particle_handler.cc.
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::insert_particle | ( | const Particle< dim, spacedim > & | particle, |
const typename Triangulation< dim, spacedim >::active_cell_iterator & | cell | ||
) |
Insert a particle into the collection of particles. Return an iterator to the new position of the particle. This function involves a copy of the particle and its properties. Note that this function is of \(O(N \log N)\) complexity for \(N\) particles.
Definition at line 256 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::insert_particles | ( | const std::multimap< typename Triangulation< dim, spacedim >::active_cell_iterator, Particle< dim, spacedim > > & | particles | ) |
Insert a number of particles into the collection of particles. This function involves a copy of the particles and their properties. Note that this function is of O(n_existing_particles + n_particles) complexity.
Definition at line 276 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::insert_particles | ( | const std::vector< Point< spacedim > > & | positions | ) |
Create and insert a number of particles into the collection of particles. This function takes a list of positions and creates a set of particles at these positions, which are then added to the local particle collection. Note that this function currently uses GridTools::compute_point_locations, which assumes all positions are within the local part of the triangulation. If one of them is not in the local domain this function will throw an exception.
Definition at line 291 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::register_additional_store_load_functions | ( | const std::function< std::size_t()> & | size_callback, |
const std::function< void *(const particle_iterator &, void *)> & | store_callback, | ||
const std::function< const void *(const particle_iterator &, const void *)> & | load_callback | ||
) |
This function allows to register three additional functions that are called every time a particle is transferred to another process (i.e. during sorting into cells, during ghost particle transfer, or during serialization of all particles).
size_callback | A function that is called when serializing particle data. The function gets no arguments and is expected to return the size of the additional data that is serialized per particle. Note that this currently implies the data size has to be the same for every particle. |
store_callback | A function that is called once per particle when serializing particle data. Arguments to the function are a particle iterator that identifies the current particle and a void pointer that points to a data block of size size_callback() in which the function can store additional data. The function is expected to return a void pointer pointing to a position right after its data block. |
load_callback | A function that is called once per particle when deserializing particle data. Arguments to the function are a particle iterator that identifies the current particle and a void pointer that points to a data block of size size_callback() in which additional data was stored by the store_callback function. The function is expected to return a void pointer pointing to a position right after its data block. |
Definition at line 869 of file particle_handler.cc.
types::particle_index Particles::ParticleHandler< dim, spacedim >::n_global_particles | ( | ) | const |
Return the total number of particles that were managed by this class the last time the update_cached_numbers() function was called. The actual number of particles may have changed since then if particles have been added or removed.
Definition at line 344 of file particle_handler.cc.
types::particle_index Particles::ParticleHandler< dim, spacedim >::n_global_max_particles_per_cell | ( | ) | const |
Return the maximum number of particles per cell the last time the update_cached_numbers() function was called.
Definition at line 353 of file particle_handler.cc.
types::particle_index Particles::ParticleHandler< dim, spacedim >::n_locally_owned_particles | ( | ) | const |
Return the number of particles in the local part of the triangulation.
Definition at line 362 of file particle_handler.cc.
types::particle_index Particles::ParticleHandler< dim, spacedim >::get_next_free_particle_index | ( | ) | const |
Return the next free particle index in the global set of particles the last time the update_cached_numbers() function was called.
Definition at line 398 of file particle_handler.cc.
unsigned int Particles::ParticleHandler< dim, spacedim >::n_properties_per_particle | ( | ) | const |
Return the number of properties each particle has.
Definition at line 371 of file particle_handler.cc.
PropertyPool & Particles::ParticleHandler< dim, spacedim >::get_property_pool | ( | ) | const |
Return a reference to the property pool that owns all particle properties, and organizes them physically.
Definition at line 407 of file particle_handler.cc.
unsigned int Particles::ParticleHandler< dim, spacedim >::n_particles_in_cell | ( | const typename Triangulation< dim, spacedim >::active_cell_iterator & | cell | ) | const |
Return the number of particles in the given cell.
Definition at line 380 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::sort_particles_into_subdomains_and_cells | ( | ) |
Find and update the cells containing each particle for all locally owned particles. If particles moved out of the local subdomain they will be sent to their new process and inserted there. After this function call every particle is either on its current process and in its current cell, or deleted (if it could not find its new process or cell).
Definition at line 446 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::exchange_ghost_particles | ( | ) |
Exchange all particles that live in cells that are ghost cells to other processes. Clears and re-populates the ghost_neighbors member variable.
Definition at line 641 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::register_store_callback_function | ( | const bool | serialization | ) |
Callback function that should be called before every refinement and when writing checkpoints. This function is used to register store_particles() with the triangulation.
Definition at line 884 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::register_load_callback_function | ( | const bool | serialization | ) |
Callback function that should be called after every refinement and after resuming from a checkpoint. This function is used to register load_particles() with the triangulation.
Definition at line 933 of file particle_handler.cc.
void Particles::ParticleHandler< dim, spacedim >::serialize | ( | Archive & | ar, |
const unsigned int | version | ||
) |
Serialize the contents of this class.
Definition at line 486 of file particle_handler.h.
|
private |
Transfer particles that have crossed subdomain boundaries to other processors. All received particles and their new cells will be appended to the received_particles
vector.
[in] | particles_to_send | All particles that should be sent and their new subdomain_ids are in this map. |
[in,out] | received_particles | Vector that stores all received particles. Note that it is not required nor checked that the list is empty, received particles are simply attached to the end of the vector. |
[in] | new_cells_for_particles | Optional vector of cell iterators with the same structure as particles_to_send . If this parameter is given it should contain the cell iterator for every particle to be send in which the particle belongs. This parameter is necessary if the cell information of the particle iterator is outdated (e.g. after particle movement). |
Definition at line 705 of file particle_handler.cc.
|
private |
Called by listener functions from Triangulation for every cell before a refinement step. All particles have to be attached to their cell to be sent around to the new processes.
Definition at line 999 of file particle_handler.cc.
|
private |
Called by listener functions after a refinement step. The local map of particles has to be read from the triangulation user_pointer.
Definition at line 1058 of file particle_handler.cc.
|
private |
Address of the triangulation to work on.
Definition at line 339 of file particle_handler.h.
|
private |
Address of the mapping to work on.
Definition at line 344 of file particle_handler.h.
|
private |
Set of particles currently living in the local domain, organized by the level/index of the cell they are in.
Definition at line 350 of file particle_handler.h.
|
private |
Set of particles that currently live in the ghost cells of the local domain, organized by the level/index of the cell they are in. These particles are equivalent to the ghost entries in distributed vectors.
Definition at line 357 of file particle_handler.h.
|
private |
This variable stores how many particles are stored globally. It is calculated by update_cached_numbers().
Definition at line 363 of file particle_handler.h.
|
private |
The maximum number of particles per cell in the global domain. This variable is important to store and load particle data during repartition and serialization of the solution. Note that the variable is only updated when it is needed, e.g. after particle movement, before/after mesh refinement, before creating a checkpoint and after resuming from a checkpoint.
Definition at line 373 of file particle_handler.h.
|
private |
This variable stores the next free particle index that is available globally in case new particles need to be generated.
Definition at line 379 of file particle_handler.h.
|
private |
This object owns and organizes the memory for all particle properties.
Definition at line 385 of file particle_handler.h.
|
private |
A function that can be registered by calling register_additional_store_load_functions. It is called when serializing particle data. The function gets no arguments and is expected to return the size of the additional data that is serialized per particle. Note that this currently implies the data size has to be the same for every particle, but it does not have to be the same for every serialization process (e.g. a serialization during particle movement might include temporary data, while a serialization after movement was finished does not need to transfer this data).
Definition at line 398 of file particle_handler.h.
|
private |
A function that can be registered by calling register_additional_store_load_functions. It is called once per particle when serializing particle data. Arguments to the function are a particle iterator that identifies the current particle and a void pointer that points to a data block of size size_callback() in which the function can store additional data. The function is expected to return a void pointer pointing to a position right after its data block.
Definition at line 411 of file particle_handler.h.
|
private |
A function that is called once per particle when deserializing particle data. Arguments to the function are a particle iterator that identifies the current particle and a void pointer that points to a data block of size size_callback() from which the function can load additional data. This block was filled by the store_callback function during serialization. This function is expected to return a void pointer pointing to a position right after its data block.
Definition at line 424 of file particle_handler.h.
|
private |
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.
Definition at line 431 of file particle_handler.h.