Reference documentation for deal.II version GIT relicensing-437-g81ec864850 2024-04-19 07:30:02+00:00
\(\newcommand{\dealvcentcolon}{\mathrel{\mathop{:}}}\) \(\newcommand{\dealcoloneq}{\dealvcentcolon\mathrel{\mkern-1.2mu}=}\) \(\newcommand{\jump}[1]{\left[\!\left[ #1 \right]\!\right]}\) \(\newcommand{\average}[1]{\left\{\!\left\{ #1 \right\}\!\right\}}\)
Loading...
Searching...
No Matches
Functions
Particles::Generators Namespace Reference

Functions

template<int dim, int spacedim = dim>
void regular_reference_locations (const Triangulation< dim, spacedim > &triangulation, const std::vector< Point< dim > > &particle_reference_locations, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
 
template<int dim, int spacedim = dim>
Particle< dim, spacedim > random_particle_in_cell (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const types::particle_index id, std::mt19937 &random_number_generator, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
 
template<int dim, int spacedim = dim>
ParticleIterator< dim, spacedim > random_particle_in_cell_insert (const typename Triangulation< dim, spacedim >::active_cell_iterator &cell, const types::particle_index id, std::mt19937 &random_number_generator, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()))
 
template<int dim, int spacedim = dim>
void probabilistic_locations (const Triangulation< dim, spacedim > &triangulation, const Function< spacedim > &probability_density_function, const bool random_cell_selection, const types::particle_index n_particles_to_create, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const unsigned int random_number_seed=5432)
 
template<int dim, int spacedim = dim>
void dof_support_points (const DoFHandler< dim, spacedim > &dof_handler, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const ComponentMask &components={}, const std::vector< std::vector< double > > &properties={})
 
template<int dim, int spacedim = dim>
void quadrature_points (const Triangulation< dim, spacedim > &triangulation, const Quadrature< dim > &quadrature, const std::vector< std::vector< BoundingBox< spacedim > > > &global_bounding_boxes, ParticleHandler< dim, spacedim > &particle_handler, const Mapping< dim, spacedim > &mapping=(ReferenceCells::get_hypercube< dim >() .template get_default_linear_mapping< dim, spacedim >()), const std::vector< std::vector< double > > &properties={})
 

Detailed Description

A namespace that contains all classes that are related to the particle generation.

Function Documentation

◆ regular_reference_locations()

template<int dim, int spacedim = dim>
void Particles::Generators::regular_reference_locations ( const Triangulation< dim, spacedim > &  triangulation,
const std::vector< Point< dim > > &  particle_reference_locations,
ParticleHandler< dim, spacedim > &  particle_handler,
const Mapping< dim, spacedim > &  mapping = (ReferenceCells::get_hypercube<dim>()           .template get_default_linear_mapping<dim, spacedim>()           ) 
)

A function that generates particles in every cell at specified particle_reference_locations. The total number of particles that is added to the particle_handler object is the number of locally owned cells of the triangulation times the number of locations in particle_reference_locations. An optional mapping argument can be used to map from particle_reference_locations to the real particle locations.

Parameters
triangulationThe triangulation associated with the particle_handler.
particle_reference_locationsA vector of positions in the unit cell. Particles will be generated in every cell at these locations.
particle_handlerThe particle handler that will take ownership of the generated particles.
mappingAn optional mapping object that is used to map reference location in the unit cell to the real cells of the triangulation. If no mapping is provided a MappingQ1 is assumed.

Definition at line 182 of file generators.cc.

◆ random_particle_in_cell()

template<int dim, int spacedim = dim>
Particle< dim, spacedim > Particles::Generators::random_particle_in_cell ( const typename Triangulation< dim, spacedim >::active_cell_iterator &  cell,
const types::particle_index  id,
std::mt19937 &  random_number_generator,
const Mapping< dim, spacedim > &  mapping = (ReferenceCells::get_hypercube<dim>()           .template get_default_linear_mapping<dim, spacedim>()           ) 
)

A function that generates one particle at a random location in cell cell and with index id. The function expects a random number generator to avoid the expensive generation and destruction of a generator for every particle and optionally takes into account a mapping for the cell. The algorithm implemented in the function is described in [88]. In short, the algorithm generates random locations within the bounding box of the cell. It then inverts the mapping to check if the generated particle is within the cell itself. This makes sure the algorithm produces statistically random locations even for nonlinear mappings and distorted cells. However, if the ratio between bounding box and cell volume becomes very large – i.e. the cells become strongly deformed, for example a pencil shaped cell that lies diagonally in the domain – then the algorithm can become very inefficient. Therefore, it only tries to find a location ni the cell a fixed number of times before throwing an error message.

Parameters
[in]cellThe cell in which a particle is generated.
[in]idThe particle index that will be assigned to the new particle.
[in,out]random_number_generatorA random number generator that will be used for the creation of th particle.
[in]mappingAn optional mapping object that is used to map reference location in the unit cell to the real cell. If no mapping is provided a MappingQ1 is assumed.

Definition at line 238 of file generators.cc.

◆ random_particle_in_cell_insert()

template<int dim, int spacedim = dim>
ParticleIterator< dim, spacedim > Particles::Generators::random_particle_in_cell_insert ( const typename Triangulation< dim, spacedim >::active_cell_iterator &  cell,
const types::particle_index  id,
std::mt19937 &  random_number_generator,
ParticleHandler< dim, spacedim > &  particle_handler,
const Mapping< dim, spacedim > &  mapping = (ReferenceCells::get_hypercube<dim>()           .template get_default_linear_mapping<dim, spacedim>()           ) 
)

A function that generates one particle at a random location in cell cell and with index id. This version of the function above immediately inserts the generated particle into the particle_handler and returns a iterator to it instead of a particle object. This avoids unnecessary copies of the particle.

Definition at line 255 of file generators.cc.

◆ probabilistic_locations()

template<int dim, int spacedim = dim>
void Particles::Generators::probabilistic_locations ( const Triangulation< dim, spacedim > &  triangulation,
const Function< spacedim > &  probability_density_function,
const bool  random_cell_selection,
const types::particle_index  n_particles_to_create,
ParticleHandler< dim, spacedim > &  particle_handler,
const Mapping< dim, spacedim > &  mapping = (ReferenceCells::get_hypercube<dim>()           .template get_default_linear_mapping<dim, spacedim>()),
const unsigned int  random_number_seed = 5432 
)

A function that generates particles randomly in the domain with a particle density according to a provided probability density function probability_density_function. The total number of particles that is added to the particle_handler object is n_particles_to_create. An optional mapping argument can be used to map from particle_reference_locations to the real particle locations. The function can compute the number of particles per cell either deterministically by computing the integral of the probability density function for each cell and creating particles accordingly (if option random_cell_selection set to false), or it can select cells randomly based on the probability density function and the cell size (if option random_cell_selection set to true). In either case the position of individual particles inside the cell is computed randomly.

The algorithm implemented in the function is described in [88].

Parameters
[in]triangulationThe triangulation associated with the particle_handler.
[in]probability_density_functionA function with non-negative values that determines the probability density of a particle to be generated in this location. The function does not need to be normalized.
[in]random_cell_selectionA bool that determines, how the number of particles per cell is computed (see the description above).
[in]n_particles_to_createThe number of particles that will be created by this function.
[in,out]particle_handlerThe particle handler that will take ownership of the generated particles.
[in]mappingAn optional mapping object that is used to map reference location in the unit cell to the real cells of the triangulation. If no mapping is provided a MappingQ1 is assumed.
[in]random_number_seedAn optional seed that determines the initial state of the random number generator. Use the same number to get a reproducible particle distribution, or a changing number (e.g. based on system time) to generate different particle distributions for each call to this function.

Definition at line 275 of file generators.cc.

◆ dof_support_points()

template<int dim, int spacedim = dim>
void Particles::Generators::dof_support_points ( const DoFHandler< dim, spacedim > &  dof_handler,
const std::vector< std::vector< BoundingBox< spacedim > > > &  global_bounding_boxes,
ParticleHandler< dim, spacedim > &  particle_handler,
const Mapping< dim, spacedim > &  mapping = (ReferenceCells::get_hypercube<dim>()           .template get_default_linear_mapping<dim, spacedim>()),
const ComponentMask components = {},
const std::vector< std::vector< double > > &  properties = {} 
)

A function that generates particles at the locations of the support points of a DoFHandler, possibly based on a different Triangulation with respect to the one used to construct the ParticleHandler. The total number of particles that is added to the particle_handler object is the number of dofs of the DoFHandler that is passed that are within the triangulation and whose components are within the ComponentMask. This function uses insert_global_particles and consequently may induce considerable mpi communication overhead.

This function is used in step-70.

Parameters
[in]dof_handlerA DOF handler that may live on another triangulation that is used to establsh the positions of the particles.
[in]global_bounding_boxesA vector that contains all the bounding boxes for all processors. This vector can be established by first using 'GridTools::compute_mesh_predicate_bounding_box()' and gathering all the bounding boxes using 'Utilities::MPI::all_gather().
[in,out]particle_handlerThe particle handler that will take ownership of the generated particles. The particles that are generated will be appended to the particles currently owned by the particle handler.
[in]mappingAn optional mapping object that is used to map the DOF locations. If no mapping is provided a MappingQ1 is assumed.
[in]componentsComponent mask that decides which subset of the support points of the dof_handler are used to generate the particles.
[in]propertiesAn optional vector of vector of properties for each particle to be inserted.

Definition at line 439 of file generators.cc.

◆ quadrature_points()

template<int dim, int spacedim = dim>
void Particles::Generators::quadrature_points ( const Triangulation< dim, spacedim > &  triangulation,
const Quadrature< dim > &  quadrature,
const std::vector< std::vector< BoundingBox< spacedim > > > &  global_bounding_boxes,
ParticleHandler< dim, spacedim > &  particle_handler,
const Mapping< dim, spacedim > &  mapping = (ReferenceCells::get_hypercube<dim>()           .template get_default_linear_mapping<dim, spacedim>()),
const std::vector< std::vector< double > > &  properties = {} 
)

A function that generates particles at the locations of the quadrature points of a Triangulation. This Triangulation can be different from the one used to construct the ParticleHandler. The total number of particles that is added to the particle_handler object is the number of cells multiplied by the number of particle reference locations which are generally constructed using a quadrature. This function uses insert_global_particles and consequently may induce considerable mpi communication overhead.

Parameters
[in]triangulationThe possibly non-matching triangulation which is used to insert the particles into the domain.
[in]quadratureA quadrature whose reference location are used to insert the particles within the cells.
[in]global_bounding_boxesA vector that contains all the bounding boxes for all processors. This vector can be established by first using 'GridTools::compute_mesh_predicate_bounding_box()' and gathering all the bounding boxes using 'Utilities::MPI::all_gather. Of course, if you are applying this function on a sequential triangulation, then this argument is simply an (outer) vector with just one element, and that one element would be the result of the function mentioned above.
[in,out]particle_handlerThe particle handler that will take ownership of the generated particles.
[in]mappingAn optional mapping object that is used to map the quadrature locations. If no mapping is provided a MappingQ1 is assumed.
[in]propertiesAn optional vector of vector of properties for each particle to be inserted.

Definition at line 473 of file generators.cc.