Reference documentation for deal.II version 9.2.0
\(\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\}}\)
Public Types | Public Member Functions | Private Member Functions | Private Attributes | List of all members
Particles::ParticleHandler< dim, spacedim > Class Template Reference

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

Public Types

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

Public Member Functions

 ParticleHandler ()
 
 ParticleHandler (const Triangulation< dim, spacedim > &tria, const Mapping< dim, spacedim > &mapping, const unsigned int n_properties=0)
 
virtual ~ParticleHandler () override=default
 
void initialize (const 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)
 
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={})
 
template<class VectorType >
std::enable_if< std::is_convertible< VectorType *, Function< spacedim > * >::value==false >::type 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<class 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_relevant_ids () const
 
unsigned int n_properties_per_particle () const
 
PropertyPoolget_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 ()
 
void register_load_callback_function (const bool serialization)
 
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 >>())
 
std::vector< char > store_particles (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status) const
 
void load_particles (const typename Triangulation< dim, spacedim >::cell_iterator &cell, const typename Triangulation< dim, spacedim >::CellStatus status, const boost::iterator_range< std::vector< char >::const_iterator > &data_range)
 

Private Attributes

SmartPointer< const 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< PropertyPoolproperty_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 handle
 
std::unique_ptr< GridTools::Cache< dim, spacedim > > triangulation_cache
 

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".

This class is used in step-70.

Definition at line 30 of file data_out.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 66 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 71 of file particle_handler.h.

Constructor & Destructor Documentation

◆ ParticleHandler() [1/2]

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

Default constructor.

Definition at line 89 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 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 106 of file particle_handler.cc.

◆ ~ParticleHandler()

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

Destructor.

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 does not clear the internal data structures, it just sets the triangulation and the mapping to be used.

Definition at line 132 of file particle_handler.cc.

◆ clear()

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

Clear all particle related data.

Definition at line 154 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 166 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 175 of file particle_handler.cc.

◆ begin() [1/2]

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

Return an iterator to the first particle.

Definition at line 229 of file particle_handler.cc.

◆ begin() [2/2]

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

Return an iterator to the first particle.

◆ end() [1/2]

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

Return an iterator past the end of the particles.

Definition at line 247 of file particle_handler.cc.

◆ end() [2/2]

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

Return an iterator past the end of the particles.

◆ begin_ghost() [1/2]

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

Return an iterator to the first ghost particle.

Definition at line 265 of file particle_handler.cc.

◆ begin_ghost() [2/2]

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

Return an iterator to the first ghost particle.

◆ end_ghost() [1/2]

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

Return an iterator past the end of the ghost particles.

Definition at line 283 of file particle_handler.cc.

◆ end_ghost() [2/2]

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

Return an iterator past the end of the ghost particles.

◆ 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.

Definition at line 313 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.

Definition at line 301 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.

Definition at line 337 of file particle_handler.cc.

◆ insert_particle()

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 347 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 371 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 391 of file particle_handler.cc.

◆ insert_global_particles()

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 = {} 
)

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 procesor. 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.

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 transfered 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.
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.
Author
Bruno Blais, Luca Heltai 2019

Definition at line 457 of file particle_handler.cc.

◆ set_particle_positions() [1/3]

template<int dim, int spacedim>
template<class VectorType >
std::enable_if< std::is_convertible< VectorType *, Function< spacedim > * >::value==false >::type Particles::ParticleHandler< dim, spacedim >::set_particle_positions ( const VectorType input_vector,
const bool  displace_particles = true 
)

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_relevant_ids(), and taking its tensor product with the index set representing the range [0, spacedim), i.e.:

IndexSet ids = particle_handler.locally_relevant_ids().
tensor_product(complete_index_set(spacedim));

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.
Authors
Luca Heltai, Bruno Blais, 2019.

Definition at line 807 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.
Authors
Bruno Blais, Luca Heltai (2019)

Definition at line 755 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 fuction 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.
Authors
Bruno Blais, Luca Heltai (2019)

Definition at line 777 of file particle_handler.cc.

◆ get_particle_positions() [1/2]

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

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.
Author
Luca Heltai, Bruno Blais, 2019.

Definition at line 830 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.
Authors
Bruno Blais, Luca Heltai (2019)

Definition at line 734 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 1395 of file particle_handler.cc.

◆ n_global_particles()

template<int dim, int spacedim>
types::particle_index Particles::ParticleHandler< dim, spacedim >::n_global_particles

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 655 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

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 664 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

Return the number of particles in the local part of the triangulation.

Definition at line 673 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

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

◆ locally_relevant_ids()

template<int dim, int spacedim>
IndexSet Particles::ParticleHandler< dim, spacedim >::locally_relevant_ids

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.
Author
Luca Heltai, Bruno Blais, 2019.

Definition at line 721 of file particle_handler.cc.

◆ n_properties_per_particle()

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

Return the number of properties each particle has.

Definition at line 682 of file particle_handler.cc.

◆ get_property_pool()

template<int dim, int spacedim>
PropertyPool & Particles::ParticleHandler< dim, spacedim >::get_property_pool

Return a reference to the property pool that owns all particle properties, and organizes them physically.

Definition at line 805 of file particle_handler.cc.

◆ n_particles_in_cell()

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

Definition at line 845 of file particle_handler.cc.

◆ exchange_ghost_particles()

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

Definition at line 1410 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 load_particles() with the triangulation. This function is used in step-70.

Definition at line 1447 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 
)

Serialize the contents of this class.

Definition at line 789 of file particle_handler.h.

◆ 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,
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>>() 
)
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.

Parameters
[in]particles_to_sendAll particles that should be sent and their new subdomain_ids are in this map.
[in,out]received_particlesVector 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_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).

Definition at line 1161 of file particle_handler.cc.

◆ store_particles()

template<int dim, int spacedim>
std::vector< char > Particles::ParticleHandler< dim, spacedim >::store_particles ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const typename Triangulation< dim, spacedim >::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 processes.

Definition at line 1512 of file particle_handler.cc.

◆ load_particles()

template<int dim, int spacedim>
void Particles::ParticleHandler< dim, spacedim >::load_particles ( const typename Triangulation< dim, spacedim >::cell_iterator &  cell,
const typename Triangulation< dim, spacedim >::CellStatus  status,
const boost::iterator_range< std::vector< char >::const_iterator > &  data_range 
)
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 1603 of file particle_handler.cc.

Member Data Documentation

◆ 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 613 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 619 of file particle_handler.h.

◆ particles

template<int dim, int spacedim = dim>
std::multimap<internal::LevelInd, Particle<dim, spacedim> > Particles::ParticleHandler< dim, spacedim >::particles
private

Set of particles currently living in the local domain, organized by the level/index of the cell they are in.

Definition at line 625 of file particle_handler.h.

◆ ghost_particles

template<int dim, int spacedim = dim>
std::multimap<internal::LevelInd, Particle<dim, spacedim> > Particles::ParticleHandler< dim, spacedim >::ghost_particles
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 632 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 638 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 648 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 654 of file particle_handler.h.

◆ property_pool

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

This object owns and organizes the memory for all particle properties.

Definition at line 660 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 673 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 685 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 698 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_store_callback_function() function and used by the register_load_callback_function() function to check where the particle data was registered in the corresponding triangulation object.

Definition at line 706 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 717 of file particle_handler.h.


The documentation for this class was generated from the following files:
IndexSet
Definition: index_set.h:74
complete_index_set
IndexSet complete_index_set(const IndexSet::size_type N)
Definition: index_set.h:1014