Reference documentation for deal.II version GIT relicensing-136-gb80d0be4af 2024-03-18 08:20: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 | Static Public Member Functions | Public Attributes | Private Types | Private Member Functions | Private Attributes | Static Private Attributes | List of all members
Particles::ParticleHandler< dim, spacedim > Class Template Reference

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

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

Classes

struct  Signals
 

Public Types

using particle_iterator = ParticleIterator< dim, spacedim >
 
using particle_iterator_range = boost::iterator_range< particle_iterator >
 
using particle_container = typename ParticleAccessor< dim, spacedim >::particle_container
 

Public Member Functions

 ParticleHandler ()
 
 ParticleHandler (const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const unsigned int n_properties=0)
 
virtual ~ParticleHandler ()
 
void initialize (const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const unsigned int n_properties=0)
 
void copy_from (const ParticleHandler< dim, spacedim > &particle_handler)
 
void clear ()
 
void clear_particles ()
 
void reserve (const std::size_t n_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 ()
 
types::particle_index n_particles_in_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell) const
 
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)
 
void remove_particles (const std::vector< particle_iterator > &particles)
 
particle_iterator insert_particle (const Particle< dim, spacedim > &particle, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
 
particle_iterator insert_particle (const Point< spacedim > &position, const Point< dim > &reference_position, const types::particle_index particle_index, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const ArrayView< const double > &properties={})
 
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)
 
std::map< unsigned int, IndexSetinsert_global_particles (const std::vector< Point< spacedim > > &positions, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, const std::vector< std::vector< double > > &properties={}, const std::vector< types::particle_index > &ids={})
 
std::map< unsigned int, IndexSetinsert_global_particles (const std::vector< Particle< dim, spacedim > > &particles, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes)
 
template<typename VectorType >
std::enable_if_t< std::is_convertible_v< VectorType *, Function< spacedim > * >==false > set_particle_positions (const VectorType &input_vector, const bool displace_particles=true)
 
void set_particle_positions (const std::vector< Point< spacedim > > &new_positions, const bool displace_particles=true)
 
void set_particle_positions (const Function< spacedim > &function, const bool displace_particles=true)
 
template<typename VectorType >
void get_particle_positions (VectorType &output_vector, const bool add_to_output_vector=false)
 
void get_particle_positions (std::vector< Point< spacedim > > &positions, const bool add_to_output_vector=false)
 
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
 
IndexSet locally_owned_particle_ids () const
 
unsigned int n_properties_per_particle () const
 
types::particle_index get_max_local_particle_index () const
 
PropertyPool< dim, spacedim > & get_property_pool () const
 
void sort_particles_into_subdomains_and_cells ()
 
void exchange_ghost_particles (const bool enable_ghost_cache=false)
 
void update_ghost_particles ()
 
void prepare_for_coarsening_and_refinement ()
 
void unpack_after_coarsening_and_refinement ()
 
void prepare_for_serialization ()
 
void deserialize ()
 
void register_store_callback_function ()
 
void register_load_callback_function (const bool serialization)
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
Subscriptor functionality

Classes derived from Subscriptor provide a facility to subscribe to this object. This is mostly used by the SmartPointer class.

void subscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
void unsubscribe (std::atomic< bool > *const validity, const std::string &identifier="") const
 
unsigned int n_subscriptions () const
 
template<typename StreamType >
void list_subscribers (StreamType &stream) const
 
void list_subscribers () const
 

Static Public Member Functions

static ::ExceptionBaseExcInUse (int arg1, std::string arg2, std::string arg3)
 
static ::ExceptionBaseExcNoSubscriber (std::string arg1, std::string arg2)
 

Public Attributes

Signals signals
 

Private Types

using map_value_type = decltype(counter_map)::value_type
 
using map_iterator = decltype(counter_map)::iterator
 

Private Member Functions

particle_iterator insert_particle (const void *&data, const typename Triangulation< dim, spacedim >::active_cell_iterator &cell)
 
particle_iterator insert_particle (const typename PropertyPool< dim, spacedim >::Handle handle, const typename Triangulation< dim, spacedim >::cell_iterator &cell)
 
void reset_particle_container (particle_container &particles)
 
void send_recv_particles (const std::map< types::subdomain_id, std::vector< particle_iterator > > &particles_to_send, 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 > >(), const bool enable_cache=false)
 
void send_recv_particles_properties_and_location (const std::map< types::subdomain_id, std::vector< particle_iterator > > &particles_to_send)
 
void connect_to_triangulation_signals ()
 
void post_mesh_change_action ()
 
void register_data_attach ()
 
void notify_ready_to_unpack (const bool serialization)
 
std::vector< char > pack_callback (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellStatus status) const
 
void unpack_callback (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range)
 
particle_container::iterator particle_container_owned_begin () const
 
particle_container::iterator particle_container_owned_end () const
 
particle_container::iterator particle_container_ghost_begin () const
 
particle_container::iterator particle_container_ghost_end () const
 
void check_no_subscribers () const noexcept
 

Private Attributes

SmartPointer< const Triangulation< dim, spacedim >, ParticleHandler< dim, spacedim > > triangulation
 
SmartPointer< const Mapping< dim, spacedim >, ParticleHandler< dim, spacedim > > mapping
 
std::unique_ptr< PropertyPool< dim, spacedim > > property_pool
 
particle_container particles
 
const particle_container::iterator owned_particles_end
 
std::vector< typename particle_container::iterator > cells_to_particle_cache
 
types::particle_index global_number_of_particles
 
types::particle_index number_of_locally_owned_particles
 
unsigned int global_max_particles_per_cell
 
types::particle_index next_free_particle_index
 
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 handle
 
const double tolerance_inside_cell = 1e-12
 
std::unique_ptr< GridTools::Cache< dim, spacedim > > triangulation_cache
 
internal::GhostParticlePartitioner< dim, spacedim > ghost_particles_cache
 
std::vector< boost::signals2::connection > tria_listeners
 
std::atomic< unsigned intcounter
 
std::map< std::string, unsigned intcounter_map
 
std::vector< std::atomic< bool > * > validity_pointers
 
const std::type_info * object_info
 

Static Private Attributes

static std::mutex mutex
 

Detailed Description

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

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". The class also includes functionality that is similar to the DoFHandler class (it knows which particles live on which cells) and the SolutionTransfer class (it knows how to transfer particles between cells and subdomains).

For examples on how to use this class to track particles, store properties on particles, and let the properties on the particles influence the finite-element solution see step-19, step-68, and step-70.

Definition at line 64 of file particle_handler.h.

Member Typedef Documentation

◆ particle_iterator

template<int dim, int spacedim = dim>
using Particles::ParticleHandler< dim, spacedim >::particle_iterator = ParticleIterator<dim, spacedim>

A type that can be used to iterate over all particles in the domain.

Definition at line 70 of file particle_handler.h.

◆ particle_iterator_range

template<int dim, int spacedim = dim>
using Particles::ParticleHandler< dim, spacedim >::particle_iterator_range = boost::iterator_range<particle_iterator>

A type that represents a range of particles.

Definition at line 75 of file particle_handler.h.

◆ particle_container

template<int dim, int spacedim = dim>
using Particles::ParticleHandler< dim, spacedim >::particle_container = typename ParticleAccessor<dim, spacedim>::particle_container

A type for the storage container for particles.

Definition at line 80 of file particle_handler.h.

◆ map_value_type

using Subscriptor::map_value_type = decltype(counter_map)::value_type
privateinherited

The data type used in counter_map.

Definition at line 229 of file subscriptor.h.

◆ map_iterator

using Subscriptor::map_iterator = decltype(counter_map)::iterator
privateinherited

The iterator type used in counter_map.

Definition at line 234 of file subscriptor.h.

Constructor & Destructor Documentation

◆ ParticleHandler() [1/2]

template<int dim, int spacedim>
Particles::ParticleHandler< dim, spacedim >::ParticleHandler ( )

Default constructor.

Definition at line 57 of file particle_handler.cc.

◆ ParticleHandler() [2/2]

template<int dim, int spacedim>
Particles::ParticleHandler< dim, spacedim >::ParticleHandler ( const 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 with 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 77 of file particle_handler.cc.

◆ ~ParticleHandler()

template<int dim, int spacedim>
Particles::ParticleHandler< dim, spacedim >::~ParticleHandler ( )
virtual

Destructor.

Definition at line 105 of file particle_handler.cc.

Member Function Documentation

◆ initialize()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::initialize ( const Triangulation< dim, spacedim > &  tria,
const Mapping< dim, spacedim > &  mapping,
const unsigned int  n_properties = 0 
)

Initialize the particle handler. This function clears the internal data structures, and sets the triangulation and the mapping to be used.

Definition at line 117 of file particle_handler.cc.

◆ copy_from()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::copy_from ( const ParticleHandler< dim, spacedim > &  particle_handler)

Copy the state of particle handler particle_handler into the current object. This will copy all particles and properties and leave this object as an identical copy of particle_handler. Existing particles in this object are deleted. Be aware that this does not copy functions that are connected to the signals of particle_handler, nor does it connect the current object's member functions to triangulation signals, which must be done by the caller if necessary, that is if the particle_handler had connected functions.

This function is expensive as it has to duplicate all data in particle_handler, and insert it into this object, which may be a significant amount of data. However, it can be useful to save the state of a particle collection at a certain point in time and reset this state later under certain conditions, for example if a timestep has to be undone and repeated.

Definition at line 148 of file particle_handler.cc.

◆ clear()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::clear ( )

Clear all particle related data but keep the handler initialized.

Definition at line 191 of file particle_handler.cc.

◆ clear_particles()

template<int dim, int spacedim>
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 204 of file particle_handler.cc.

◆ reserve()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::reserve ( const std::size_t  n_particles)

This function can be used to preemptively reserve memory for particle data. Calling this function before inserting particles will reduce memory allocations and therefore increase the performance. Calling this function is optional; if memory is not already allocated it will be allocated automatically during the insertion. It is recommended to use this function if you know the number of particles that will be inserted, but cannot use one of the collective particle insertion functions.

Parameters
n_particlesNumber of particles to reserve memory for. Note that this is the total number of particles to be stored, not the number of particles to be newly inserted.

Definition at line 226 of file particle_handler.cc.

◆ update_cached_numbers()

template<int dim, int spacedim>
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 261 of file particle_handler.cc.

◆ begin() [1/2]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::begin ( ) const
inline

Return an iterator to the first locally owned particle.

Definition at line 1199 of file particle_handler.h.

◆ begin() [2/2]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::begin ( )
inline

Return an iterator to the first locally owned particle.

Definition at line 1208 of file particle_handler.h.

◆ end() [1/2]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::end ( ) const
inline

Return an iterator past the end of the particles.

Definition at line 1219 of file particle_handler.h.

◆ end() [2/2]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::end ( )
inline

Return an iterator past the end of the particles.

Definition at line 1228 of file particle_handler.h.

◆ begin_ghost() [1/2]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::begin_ghost ( ) const
inline

Return an iterator to the first ghost particle.

Definition at line 1237 of file particle_handler.h.

◆ begin_ghost() [2/2]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::begin_ghost ( )
inline

Return an iterator to the first ghost particle.

Definition at line 1246 of file particle_handler.h.

◆ end_ghost() [1/2]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::end_ghost ( ) const
inline

Return an iterator past the end of the ghost particles.

Definition at line 1257 of file particle_handler.h.

◆ end_ghost() [2/2]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::end_ghost ( )
inline

Return an iterator past the end of the ghost particles.

Definition at line 1266 of file particle_handler.h.

◆ n_particles_in_cell()

template<int dim, int spacedim>
types::particle_index Particles::ParticleHandler< dim, spacedim >::n_particles_in_cell ( const typename Triangulation< dim, spacedim >::active_cell_iterator &  cell) const

Return the number of particles that live on the given cell.

Definition at line 409 of file particle_handler.cc.

◆ particles_in_cell() [1/2]

template<int dim, int spacedim>
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.

The number of elements in the returned range equals what the n_particles_in_cell() function returns.

Definition at line 449 of file particle_handler.cc.

◆ particles_in_cell() [2/2]

template<int dim, int spacedim>
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.

The number of elements in the returned range equals what the n_particles_in_cell() function returns.

Definition at line 437 of file particle_handler.cc.

◆ remove_particle()

template<int dim, int spacedim = dim>
void Particles::ParticleHandler< dim, spacedim >::remove_particle ( const particle_iterator particle)

Remove a particle pointed to by the iterator. Note that particle and all iterators that point to other particles in the same cell as particle will be invalidated during this call.

Definition at line 488 of file particle_handler.cc.

◆ remove_particles()

template<int dim, int spacedim = dim>
void Particles::ParticleHandler< dim, spacedim >::remove_particles ( const std::vector< particle_iterator > &  particles)

Remove a vector of particles indicated by the particle iterators. The iterators and all other particle iterators are invalidated during the function call.

Definition at line 525 of file particle_handler.cc.

◆ insert_particle() [1/4]

template<int dim, int spacedim>
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 578 of file particle_handler.cc.

◆ insert_particle() [2/4]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::insert_particle ( const Point< spacedim > &  position,
const Point< dim > &  reference_position,
const types::particle_index  particle_index,
const typename Triangulation< dim, spacedim >::active_cell_iterator &  cell,
const ArrayView< const double > &  properties = {} 
)

Insert a particle into the collection of particles given all the properties necessary for creating a particle. This function is used to efficiently generate particles without the detour through a Particle object.

Parameters
[in]positionInitial position of the particle in real space.
[in]reference_positionInitial position of the particle in the coordinate system of the reference cell.
[in]particle_indexGlobally unique identifier for this particle.
[in]cellThe cell in which the particle is located.
[in]propertiesAn optional ArrayView that describes the particle properties. If given this has to be of size n_properties_per_particle().

Definition at line 649 of file particle_handler.cc.

◆ insert_particles() [1/2]

template<int dim, int spacedim>
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 683 of file particle_handler.cc.

◆ insert_particles() [2/2]

template<int dim, int spacedim>
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 699 of file particle_handler.cc.

◆ insert_global_particles() [1/2]

template<int dim, int spacedim>
std::map< unsigned int, IndexSet > Particles::ParticleHandler< dim, spacedim >::insert_global_particles ( const std::vector< Point< spacedim > > &  positions,
const std::vector< std::vector< BoundingBox< spacedim > > > &  global_bounding_boxes,
const std::vector< std::vector< double > > &  properties = {},
const std::vector< types::particle_index > &  ids = {} 
)

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 distributed and added to the local particle collection of a processor. Note that this function uses GridTools::distributed_compute_point_locations(). Consequently, it can require intense communications between the processors. This function is used in step-70.

This function figures out what mpi process owns the points that do not fall within the locally owned part of the triangulation, it sends to that process the points passed to this function on this process, and receives the points that fall within the locally owned cells of the triangulation from whoever received them as input.

In order to keep track of what mpi process received what points, a map from mpi process to IndexSet is returned by the function. This IndexSet contains the local indices of the points that were passed to this function on the calling mpi process, and that falls within the part of triangulation owned by this mpi process.

The ids of the resulting particles are assigned from the optional argument ids. If the vector of ids is empty, then the ids are computed automatically from the get_next_free_particle_index() onward. For example, if the method get_next_free_particle_index() returns n0, calling this function with two MPI processes each adding n1 and n2 particles will result in the n1 particles added by process zero having ids equal to [n0,n0+n1), and the n2 particles added by process one having ids [n0+n1, n0+n1+n2).

Parameters
[in]positionsA vector of points that do not need to be on the local processor, but have to be in the triangulation that is associated with this ParticleHandler object.
[in]global_bounding_boxesA vector of vectors of bounding boxes. The bounding boxes global_bboxes[rk] describe which part of the mesh is locally owned by the mpi process with rank rk. The local description can be obtained from GridTools::compute_mesh_predicate_bounding_box(), and the global one can be obtained by passing the local ones to Utilities::MPI::all_gather().
[in]properties(Optional) A vector of vector of properties associated with each local point. The size of the vector should be either zero (no properties will be transferred nor attached to the generated particles) or it should be a vector of positions.size() vectors of size n_properties_per_particle(). Notice that this function call will transfer the properties from the local mpi process to the final mpi process that will own each of the particles, and it may therefore be communication intensive.
[in]ids(Optional) A vector of ids to associate to each particle. If the vector is empty, the ids are assigned as a continuous range from the first available index, as documented above. If the vector is not empty, then its size must match the size of the positions vector.
Returns
A map from owner to IndexSet, that contains the local indices of the points that were passed to this function on the calling mpi process, and that falls within the part of triangulation owned by this mpi process.

Definition at line 763 of file particle_handler.cc.

◆ insert_global_particles() [2/2]

template<int dim, int spacedim>
std::map< unsigned int, IndexSet > Particles::ParticleHandler< dim, spacedim >::insert_global_particles ( const std::vector< Particle< dim, spacedim > > &  particles,
const std::vector< std::vector< BoundingBox< spacedim > > > &  global_bounding_boxes 
)

Insert a number of particles into the collection of particles. This function takes a list of particles for which we don't know the associated cell iterator, and distributes them to the correct local particle collection of a processor, by unpacking the locations, figuring out where to send the particles by calling GridTools::distributed_compute_point_locations(), and sending the particles to the corresponding process.

In order to keep track of what mpi process received what particles, a map from mpi process to IndexSet is returned by the function. This IndexSet contains the local indices of the particles that were passed to this function on the calling mpi process, and that falls within the part of the triangulation owned by this mpi process.

Parameters
[in]particlesA vector of particles that do not need to be on the local processor.
[in]global_bounding_boxesA vector of vectors of bounding boxes. The bounding boxes global_bboxes[rk] describe which part of the mesh is locally owned by the mpi process with rank rk. The local description can be obtained from GridTools::compute_mesh_predicate_bounding_box(), and the global one can be obtained by passing the local ones to Utilities::MPI::all_gather().
Returns
A map from owner to IndexSet, that contains the local indices of the points that were passed to this function on the calling mpi process, and that falls within the part of triangulation owned by this mpi process.

Definition at line 991 of file particle_handler.cc.

◆ set_particle_positions() [1/3]

template<int dim, int spacedim>
template<typename VectorType >
std::enable_if_t< std::is_convertible_v< VectorType *, Function< spacedim > * >==false > Particles::ParticleHandler< dim, spacedim >::set_particle_positions ( const VectorType &  input_vector,
const bool  displace_particles = true 
)
inline

Set the position of the particles by using the values contained in the vector input_vector.

Template Parameters
VectorTypeAny of the parallel distributed vectors supported by the library.

The vector input_vector should have read access to the indices created by extracting the locally relevant ids with locally_owned_particle_ids(), and taking its tensor product with the index set representing the range [0, spacedim), i.e.:

IndexSet ids = particle_handler.locally_owned_particle_ids().
tensor_product(complete_index_set(spacedim));
IndexSet complete_index_set(const IndexSet::size_type N)
Definition index_set.h:1193

The position of the particle with global index id is read from spacedim consecutive entries starting from input_vector[id*spacedim].

Notice that it is not necessary that the input_vector owns those indices, however it has to have read access to them (i.e., it can be a distributed vector with ghost entries).

If the argument displace_particles is set to false, then the new position taken from the values contained in input_vector, replacing the previously stored particle position. By default, the particles are displaced by the amount contained in the input_vector, i.e., the contents of the vector are considered offsets that are added to the previous position.

After setting the new position, this function calls internally the method sort_particles_into_subdomains_and_cells(). You should make sure you satisfy the requirements of that function.

Parameters
[in]input_vectorA parallel distributed vector containing the displacement to apply to each particle, or their new absolute position.
[in]displace_particlesControl if the input_vector should be interpreted as a displacement vector, or a vector of absolute positions.

Definition at line 1340 of file particle_handler.h.

◆ set_particle_positions() [2/3]

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::set_particle_positions ( const std::vector< Point< spacedim > > &  new_positions,
const bool  displace_particles = true 
)

Set the position of the particles within the particle handler using a vector of points. The new set of point defined by the vector has to be sufficiently close to the original one to ensure that the sort_particles_into_subdomains_and_cells() function manages to find the new cells in which the particles belong.

Points are numbered in the same way they are traversed locally by the ParticleHandler. A typical way to use this method, is to first call the get_particle_positions() function, and then modify the resulting vector.

Parameters
[in]new_positionsA vector of points of dimension particle_handler.n_locally_owned_particles()
[in]displace_particlesWhen true, this function adds the value of the vector of points to the current position of the particle, thus displacing them by the amount given by the function. When false, the position of the particle is replaced by the value in the vector.

Definition at line 1119 of file particle_handler.cc.

◆ set_particle_positions() [3/3]

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::set_particle_positions ( const Function< spacedim > &  function,
const bool  displace_particles = true 
)

Set the position of the particles within the particle handler using a function with spacedim components. The new set of point defined by the function has to be sufficiently close to the original one to ensure that the sort_particles_into_subdomains_and_cells algorithm manages to find the new cells in which the particles belong.

The function is evaluated at the current location of the particles.

Parameters
[in]functionA function that has n_components==spacedim that describes either the displacement or the new position of the particles as a function of the current location of the particle.
[in]displace_particlesWhen true, this function adds the results of the function to the current position of the particle, thus displacing them by the amount given by the function. When false, the position of the particle is replaced by the value of the function.

Definition at line 1142 of file particle_handler.cc.

◆ get_particle_positions() [1/2]

template<int dim, int spacedim>
template<typename VectorType >
void Particles::ParticleHandler< dim, spacedim >::get_particle_positions ( VectorType &  output_vector,
const bool  add_to_output_vector = false 
)
inline

Read the position of the particles and store them into the distributed vector output_vector. By default the output_vector is overwritten by this operation, but you can add to its entries by setting add_to_output_vector to true.

Template Parameters
VectorTypeAny of the parallel distributed vectors supported by the library.

This is the reverse operation of the set_particle_positions() function. The position of the particle with global index id is written to spacedim consecutive entries starting from output_vector[id*spacedim].

Notice that, if you use a distributed vector type, it is not necessary for the output_vector to own the entries corresponding to the indices that will be written. However you should keep in mind that this requires a global communication to distribute the entries above to their respective owners.

Parameters
[in,out]output_vectorA parallel distributed vector containing the positions of the particles, or updated with the positions of the particles.
[in]add_to_output_vectorControl if the function should set the entries of the output_vector or if should add to them.

Definition at line 1366 of file particle_handler.h.

◆ get_particle_positions() [2/2]

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::get_particle_positions ( std::vector< Point< spacedim > > &  positions,
const bool  add_to_output_vector = false 
)

Gather the position of the particles within the particle handler in a vector of points. The order of the points is the same on would obtain by iterating over all (local) particles, and querying their locations.

Parameters
[in,out]positionsA vector preallocated at size particle_handler.n_locally_owned_articles and whose points will become the positions of the locally owned particles
[in]add_to_output_vectorWhen true, the value of the point of the particles is added to the positions vector. When false, the value of the points in the positions vector are replaced by the position of the particles.

Definition at line 1098 of file particle_handler.cc.

◆ register_additional_store_load_functions()

template<int dim, int spacedim>
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).

Parameters
size_callbackA 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_callbackA 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_callbackA 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 2069 of file particle_handler.cc.

◆ n_global_particles()

template<int dim, int spacedim>
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.

Returns
Total number of particles in simulation.

Definition at line 1028 of file particle_handler.cc.

◆ n_global_max_particles_per_cell()

template<int dim, int spacedim>
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.

Returns
Maximum number of particles in one cell in simulation.

Definition at line 1037 of file particle_handler.cc.

◆ n_locally_owned_particles()

template<int dim, int spacedim>
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 1046 of file particle_handler.cc.

◆ get_next_free_particle_index()

template<int dim, int spacedim>
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 1064 of file particle_handler.cc.

◆ locally_owned_particle_ids()

template<int dim, int spacedim>
IndexSet Particles::ParticleHandler< dim, spacedim >::locally_owned_particle_ids ( ) const

Extract an IndexSet with global dimensions equal to get_next_free_particle_index(), containing the locally owned particle indices.

This function can be used to construct distributed vectors and matrices to manipulate particles using linear algebra operations.

Notice that it is the user's responsibility to guarantee that particle indices are unique, and no check is performed to verify that this is the case, nor that the union of all IndexSet objects on each mpi process is complete.

Returns
An IndexSet of size get_next_free_particle_index(), containing n_locally_owned_particle() indices.

Definition at line 1073 of file particle_handler.cc.

◆ n_properties_per_particle()

template<int dim, int spacedim>
unsigned int Particles::ParticleHandler< dim, spacedim >::n_properties_per_particle ( ) const

Return the number of properties each particle has.

Definition at line 1055 of file particle_handler.cc.

◆ get_max_local_particle_index()

template<int dim, int spacedim>
types::particle_index Particles::ParticleHandler< dim, spacedim >::get_max_local_particle_index ( ) const

Return one past the largest local index (in MPI-local index space) returned by ParticleAccessor::get_local_index(). This number can be larger than locally_owned_particle_ids().n_elements(), because local indices are not necessarily forming a contiguous range and can contain gaps in the range 0 to get_max_local_particle_index(). As a consequence, the number is not updated upon calls to remove_particle() or similar functions, and refreshed only as new particles get added or in sort_particles_into_subdomains_and_cells().

This function is appropriate for resizing vectors working with ParticleAccessor::get_local_index().

Definition at line 1089 of file particle_handler.cc.

◆ get_property_pool()

template<int dim, int spacedim>
PropertyPool< dim, spacedim > & 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 1169 of file particle_handler.cc.

◆ sort_particles_into_subdomains_and_cells()

template<int dim, int spacedim>
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).

The user may attach a function to the signal Particles::ParticleHandler::Signals::particle_lost(). The signal is triggered whenever a particle is deleted, and the connected functions are called passing an iterator to the particle in question, and its last known cell association.

Definition at line 1209 of file particle_handler.cc.

◆ exchange_ghost_particles()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::exchange_ghost_particles ( const bool  enable_ghost_cache = false)

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 1518 of file particle_handler.cc.

◆ update_ghost_particles()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::update_ghost_particles ( )

Update all particles that live in cells that are ghost cells to other processes. In this context, update means to update the location and the properties of the ghost particles assuming that the ghost particles have not changed cells. Consequently, this will not update the reference location of the particles.

Definition at line 1627 of file particle_handler.cc.

◆ prepare_for_coarsening_and_refinement()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::prepare_for_coarsening_and_refinement ( )

This function prepares the particle handler for a coarsening and refinement cycle, by storing the necessary information to transfer particles to their new cells. The implementation depends on the triangulation type that is connected to the particle handler and differs between shared and distributed triangulations. This function should be used like the corresponding function with the same name in the SolutionTransfer() class.

Definition at line 2153 of file particle_handler.cc.

◆ unpack_after_coarsening_and_refinement()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::unpack_after_coarsening_and_refinement ( )

This function unpacks the particle data after a coarsening and refinement cycle, by reading the necessary information to transfer particles to their new cells. This function should be used like the SolutionTransfer::interpolate() function after mesh refinement has finished. Note that this function requires a working mapping, i.e. if you use a mapping class that requires setup after a mesh refinement (e.g. MappingQCache(), or MappingEulerian()), the mapping has to be ready before you can call this function.

Note
It is important to note, that if you do not call this function after a refinement/coarsening operation (and a repartitioning operation in a distributed triangulation) and the particle handler contained particles before, the particle handler will not be usable any more. Not calling this function after refinement would be similar to trying to access a solution vector after mesh refinement without first using a SolutionTransfer class to transfer the vector to the new mesh.

Definition at line 2212 of file particle_handler.cc.

◆ prepare_for_serialization()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::prepare_for_serialization ( )

This function prepares the particle handler for serialization. This needs to be done before calling Triangulation::save(). This function should be used like the corresponding function with the same name in the SolutionTransfer() class.

Note
It is important to note, that if you do not call this function before a serialization operation in a distributed triangulation, the particles will not be serialized, and therefore will disappear after deserialization.

Definition at line 2162 of file particle_handler.cc.

◆ deserialize()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::deserialize ( )

Execute the deserialization of the particle data. This needs to be done after calling Triangulation::load(). The data must have been stored before the serialization of the triangulation using the prepare_for_serialization() function. This function should be used like the corresponding function with the same name in the SolutionTransfer class.

Definition at line 2222 of file particle_handler.cc.

◆ register_store_callback_function()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::register_store_callback_function ( )

Callback function that should be called before every refinement and when writing checkpoints. This function is used to register pack_callback() with the triangulation. This function is used in step-70.

Deprecated:
Please use prepare_for_coarsening_and_refinement() or prepare_for_serialization() instead. See there for further information about the purpose of this function.

Definition at line 2171 of file particle_handler.cc.

◆ register_load_callback_function()

template<int dim, int spacedim>
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 unpack_callback() with the triangulation. This function is used in step-70.

Deprecated:
Please use unpack_after_coarsening_and_refinement() or deserialize() instead. See there for further information about the purpose of this function.

Definition at line 2231 of file particle_handler.cc.

◆ serialize()

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

Serialize the contents of this class using the BOOST serialization library.

Definition at line 1323 of file particle_handler.h.

◆ insert_particle() [3/4]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::insert_particle ( const void *&  data,
const typename Triangulation< dim, spacedim >::active_cell_iterator &  cell 
)
private

Insert a particle into the collection of particles from a raw data pointer. This function is used for shipping particles during serialization and refinement and not intended for use outside of this class. Return an iterator to the new position of the particle.

Definition at line 624 of file particle_handler.cc.

◆ insert_particle() [4/4]

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_iterator Particles::ParticleHandler< dim, spacedim >::insert_particle ( const typename PropertyPool< dim, spacedim >::Handle  handle,
const typename Triangulation< dim, spacedim >::cell_iterator &  cell 
)
private

Perform the local insertion operation into the particle container. This function is used in the higher-level functions inserting particles.

Definition at line 593 of file particle_handler.cc.

◆ reset_particle_container()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::reset_particle_container ( particle_container particles)
private

Delete all entries in the particles container, and then set the three anchor entries used to distinguish between owned and ghost cells: We add one item to the front, one between and one as the last element in the particle container to be able to iterate across particles without if statements, solely relying on ParticleAccessor::operator== to terminate operations, and using the cell_iterator inside the particle container to check for valid states.

Definition at line 235 of file particle_handler.cc.

◆ send_recv_particles()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::send_recv_particles ( const std::map< types::subdomain_id, std::vector< particle_iterator > > &  particles_to_send,
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>>(),
const bool  enable_cache = false 
)
private

Transfer particles that have crossed subdomain boundaries to other processors. All received particles and their new cells will be appended to the class variable particles at the right slot.

Parameters
[in]particles_to_sendAll particles that should be sent and their new subdomain_ids are in this map.
[in]new_cells_for_particlesOptional 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).
[in]enable_cacheOptional bool that enables updating the ghost particles without rebuilding them from scratch by building a cache of type GhostParticlePartitioner, which stores the necessary information to update the ghost particles. Once this cache is built, the ghost particles can be updated by a call to send_recv_particles_properties_and_location().

Definition at line 1660 of file particle_handler.cc.

◆ send_recv_particles_properties_and_location()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::send_recv_particles_properties_and_location ( const std::map< types::subdomain_id, std::vector< particle_iterator > > &  particles_to_send)
private

Transfer ghost particles' position and properties assuming that the particles have not changed cells. This routine uses the GhostParticlePartitioner as a caching structure to know which particles are ghost to other processes, and where they need to be sent. It inherently assumes that particles cannot have changed cell, and writes the result back to the particles member variable.

Parameters
[in]particles_to_sendAll particles for which information should be sent and their new subdomain_ids are in this map.

Definition at line 1955 of file particle_handler.cc.

◆ connect_to_triangulation_signals()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::connect_to_triangulation_signals ( )
private

Connect the particle handler to the relevant triangulation signals to appropriately react to changes in the underlying triangulation.

Definition at line 2083 of file particle_handler.cc.

◆ post_mesh_change_action()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::post_mesh_change_action ( )
private

Function that gets called by the triangulation signals if the structure of the mesh has changed. This function is responsible for resizing the particle container if no particles are stored, i.e. if the usual call to prepare_for_..., and transfer_particles_after_... functions is not necessary.

Definition at line 2125 of file particle_handler.cc.

◆ register_data_attach()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::register_data_attach ( )
private

Function that should be called before every refinement and when writing checkpoints. This function is used to register pack_callback() with the triangulation.

Definition at line 2180 of file particle_handler.cc.

◆ notify_ready_to_unpack()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::notify_ready_to_unpack ( const bool  serialization)
private

Callback function that should be called after every refinement and after resuming from a checkpoint. This function is used to call the notify_ready_to_unpack() function of the triangulation and hand over the unpack_callback() function, which will unpack the particle data for each cell.

Definition at line 2241 of file particle_handler.cc.

◆ pack_callback()

template<int dim, int spacedim>
std::vector< char > Particles::ParticleHandler< dim, spacedim >::pack_callback ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const CellStatus  status 
) const
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 cell and owning process.

Definition at line 2294 of file particle_handler.cc.

◆ unpack_callback()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::unpack_callback ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const CellStatus  status,
const boost::iterator_range< std::vector< char >::const_iterator > &  data_range 
)
private

Called by listener functions after a refinement step for each cell to unpack the particle data and transfer it to the local container.

Definition at line 2350 of file particle_handler.cc.

◆ particle_container_owned_begin()

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_container::iterator Particles::ParticleHandler< dim, spacedim >::particle_container_owned_begin ( ) const
inlineprivate

Internal function returning an iterator to the begin of the container for owned particles.

Definition at line 1275 of file particle_handler.h.

◆ particle_container_owned_end()

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_container::iterator Particles::ParticleHandler< dim, spacedim >::particle_container_owned_end ( ) const
inlineprivate

Internal function returning an iterator to the end of the container for owned particles.

Definition at line 1289 of file particle_handler.h.

◆ particle_container_ghost_begin()

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_container::iterator Particles::ParticleHandler< dim, spacedim >::particle_container_ghost_begin ( ) const
inlineprivate

Internal function returning an iterator to the begin of the container for ghost particles.

Definition at line 1299 of file particle_handler.h.

◆ particle_container_ghost_end()

template<int dim, int spacedim>
ParticleHandler< dim, spacedim >::particle_container::iterator Particles::ParticleHandler< dim, spacedim >::particle_container_ghost_end ( ) const
inlineprivate

Internal function returning an iterator to the end of the container for ghost particles.

Definition at line 1310 of file particle_handler.h.

◆ subscribe()

void Subscriptor::subscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Subscribes a user of the object by storing the pointer validity. The subscriber may be identified by text supplied as identifier.

Definition at line 135 of file subscriptor.cc.

◆ unsubscribe()

void Subscriptor::unsubscribe ( std::atomic< bool > *const  validity,
const std::string &  identifier = "" 
) const
inherited

Unsubscribes a user from the object.

Note
The identifier and the validity pointer must be the same as the one supplied to subscribe().

Definition at line 155 of file subscriptor.cc.

◆ n_subscriptions()

unsigned int Subscriptor::n_subscriptions ( ) const
inlineinherited

Return the present number of subscriptions to this object. This allows to use this class for reference counted lifetime determination where the last one to unsubscribe also deletes the object.

Definition at line 300 of file subscriptor.h.

◆ list_subscribers() [1/2]

template<typename StreamType >
void Subscriptor::list_subscribers ( StreamType &  stream) const
inlineinherited

List the subscribers to the input stream.

Definition at line 317 of file subscriptor.h.

◆ list_subscribers() [2/2]

void Subscriptor::list_subscribers ( ) const
inherited

List the subscribers to deallog.

Definition at line 203 of file subscriptor.cc.

◆ check_no_subscribers()

void Subscriptor::check_no_subscribers ( ) const
privatenoexceptinherited

Check that there are no objects subscribing to this object. If this check passes then it is safe to destroy the current object. It this check fails then this function will either abort or print an error message to deallog (by using the AssertNothrow mechanism), but will not throw an exception.

Note
Since this function is just a consistency check it does nothing in release mode.
If this function is called when there is an uncaught exception then, rather than aborting, this function prints an error message to the standard error stream and returns.

Definition at line 52 of file subscriptor.cc.

Member Data Documentation

◆ signals

template<int dim, int spacedim = dim>
Signals Particles::ParticleHandler< dim, spacedim >::signals
mutable

Signals for the events that a particle handler can notify the calling application about.

Definition at line 869 of file particle_handler.h.

◆ triangulation

template<int dim, int spacedim = dim>
SmartPointer<const Triangulation<dim, spacedim>, ParticleHandler<dim, spacedim> > Particles::ParticleHandler< dim, spacedim >::triangulation
private

Address of the triangulation to work on.

Definition at line 909 of file particle_handler.h.

◆ mapping

template<int dim, int spacedim = dim>
SmartPointer<const Mapping<dim, spacedim>, ParticleHandler<dim, spacedim> > Particles::ParticleHandler< dim, spacedim >::mapping
private

Address of the mapping to work on.

Definition at line 915 of file particle_handler.h.

◆ property_pool

template<int dim, int spacedim = dim>
std::unique_ptr<PropertyPool<dim, spacedim> > Particles::ParticleHandler< dim, spacedim >::property_pool
private

This object owns and organizes the memory for all particle properties. Since particles reference the property pool, the latter has to be destroyed after the particles are destroyed. This is achieved by making sure the property_pool member variable precedes the declaration of the particles member variable.

Definition at line 925 of file particle_handler.h.

◆ particles

template<int dim, int spacedim = dim>
particle_container Particles::ParticleHandler< dim, spacedim >::particles
private

Set of particles currently living in the locally owned or ghost cells.

Definition at line 930 of file particle_handler.h.

◆ owned_particles_end

template<int dim, int spacedim = dim>
const particle_container::iterator Particles::ParticleHandler< dim, spacedim >::owned_particles_end
private

Iterator to the end of the list elements of particle_container which belong to locally owned elements. Made const to avoid accidental modification.

Definition at line 937 of file particle_handler.h.

◆ cells_to_particle_cache

template<int dim, int spacedim = dim>
std::vector<typename particle_container::iterator> Particles::ParticleHandler< dim, spacedim >::cells_to_particle_cache
private

List from the active cells on the present MPI process to positions in either owned_particles or ghost_particles for fast \(\mathcal O(1)\) access to the particles of a cell.

Definition at line 944 of file particle_handler.h.

◆ global_number_of_particles

template<int dim, int spacedim = dim>
types::particle_index Particles::ParticleHandler< dim, spacedim >::global_number_of_particles
private

This variable stores how many particles are stored globally. It is calculated by update_cached_numbers().

Definition at line 950 of file particle_handler.h.

◆ number_of_locally_owned_particles

template<int dim, int spacedim = dim>
types::particle_index Particles::ParticleHandler< dim, spacedim >::number_of_locally_owned_particles
private

This variable stores how many particles are locally owned. It is calculated by update_cached_numbers().

Definition at line 956 of file particle_handler.h.

◆ global_max_particles_per_cell

template<int dim, int spacedim = dim>
unsigned int Particles::ParticleHandler< dim, spacedim >::global_max_particles_per_cell
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 966 of file particle_handler.h.

◆ next_free_particle_index

template<int dim, int spacedim = dim>
types::particle_index Particles::ParticleHandler< dim, spacedim >::next_free_particle_index
private

This variable stores the next free particle index that is available globally in case new particles need to be generated.

Definition at line 972 of file particle_handler.h.

◆ size_callback

template<int dim, int spacedim = dim>
std::function<std::size_t()> Particles::ParticleHandler< dim, spacedim >::size_callback
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 985 of file particle_handler.h.

◆ store_callback

template<int dim, int spacedim = dim>
std::function<void *(const particle_iterator &, void *)> Particles::ParticleHandler< dim, spacedim >::store_callback
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 997 of file particle_handler.h.

◆ load_callback

template<int dim, int spacedim = dim>
std::function<const void *(const particle_iterator &, const void *)> Particles::ParticleHandler< dim, spacedim >::load_callback
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 1010 of file particle_handler.h.

◆ handle

template<int dim, int spacedim = dim>
unsigned int Particles::ParticleHandler< dim, spacedim >::handle
private

This variable is set by the register_data_attach() function and used by the notify_ready_to_unpack() function to check where the particle data was registered in the corresponding triangulation object.

Definition at line 1018 of file particle_handler.h.

◆ tolerance_inside_cell

template<int dim, int spacedim = dim>
const double Particles::ParticleHandler< dim, spacedim >::tolerance_inside_cell = 1e-12
private

Tolerance to be used for GeometryInfo<dim>::is_inside_cell().

Definition at line 1023 of file particle_handler.h.

◆ triangulation_cache

template<int dim, int spacedim = dim>
std::unique_ptr<GridTools::Cache<dim, spacedim> > Particles::ParticleHandler< dim, spacedim >::triangulation_cache
private

The GridTools::Cache is used to store the information about the vertex_to_cells set and the vertex_to_cell_centers vectors to prevent recomputing them every time we sort_into_subdomain_and_cells(). This cache is automatically updated when the triangulation has changed. This cache is stored within a unique pointer because the particle handler has a constructor that enables it to be constructed without a triangulation. The cache does not have such a constructor.

Definition at line 1034 of file particle_handler.h.

◆ ghost_particles_cache

template<int dim, int spacedim = dim>
internal::GhostParticlePartitioner<dim, spacedim> Particles::ParticleHandler< dim, spacedim >::ghost_particles_cache
private

Cache structure used to store the elements which are required to exchange the particle information (location and properties) across processors in order to update the ghost particles. This structure is only used to update the ghost particles.

Definition at line 1098 of file particle_handler.h.

◆ tria_listeners

template<int dim, int spacedim = dim>
std::vector<boost::signals2::connection> Particles::ParticleHandler< dim, spacedim >::tria_listeners
private

A list of connections with which this object connects to the triangulation to get information about when the triangulation changes.

Definition at line 1111 of file particle_handler.h.

◆ counter

std::atomic<unsigned int> Subscriptor::counter
mutableprivateinherited

Store the number of objects which subscribed to this object. Initially, this number is zero, and upon destruction it shall be zero again (i.e. all objects which subscribed should have unsubscribed again).

The creator (and owner) of an object is counted in the map below if HE manages to supply identification.

We use the mutable keyword in order to allow subscription to constant objects also.

This counter may be read from and written to concurrently in multithreaded code: hence we use the std::atomic class template.

Definition at line 218 of file subscriptor.h.

◆ counter_map

std::map<std::string, unsigned int> Subscriptor::counter_map
mutableprivateinherited

In this map, we count subscriptions for each different identification string supplied to subscribe().

Definition at line 224 of file subscriptor.h.

◆ validity_pointers

std::vector<std::atomic<bool> *> Subscriptor::validity_pointers
mutableprivateinherited

In this vector, we store pointers to the validity bool in the SmartPointer objects that subscribe to this class.

Definition at line 240 of file subscriptor.h.

◆ object_info

const std::type_info* Subscriptor::object_info
mutableprivateinherited

Pointer to the typeinfo object of this object, from which we can later deduce the class name. Since this information on the derived class is neither available in the destructor, nor in the constructor, we obtain it in between and store it here.

Definition at line 248 of file subscriptor.h.

◆ mutex

std::mutex Subscriptor::mutex
staticprivateinherited

A mutex used to ensure data consistency when accessing the mutable members of this class. This lock is used in the subscribe() and unsubscribe() functions, as well as in list_subscribers().

Definition at line 271 of file subscriptor.h.


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