ASPECT
Classes | Public Member Functions | Static Public Member Functions | Private Member Functions | Private Attributes | List of all members
aspect::Particle::World< dim > Class Template Reference
Inheritance diagram for aspect::Particle::World< dim >:
Inheritance graph
[legend]

Classes

struct  ParticleLoadBalancing
 

Public Member Functions

 World ()
 
 ~World () override
 
void initialize ()
 
const Property::Manager< dim > & get_property_manager () const
 
const Particles::ParticleHandler< dim > & get_particle_handler () const
 
Particles::ParticleHandler< dim > & get_particle_handler ()
 
void copy_particle_handler (const Particles::ParticleHandler< dim > &from_particle_handler, Particles::ParticleHandler< dim > &to_particle_handler) const
 
void backup_particles ()
 Stores a copy of the particle handler in particle_handler_backup. This copy can be used to restore the position and properties of the particles for example after advection solver iterations of the iterated advection scheme or when a timestep has to be repeated. More...
 
void restore_particles ()
 Restores the particle handler particle_handler based on the copy in particle_handler_backup. This restores the position and properties of the particles to those in the copy and can be used for example after advection solver iterations of the iterated advection scheme or when a timestep has to be repeated. More...
 
void setup_initial_state ()
 
const Interpolator::Interface< dim > & get_interpolator () const
 
void generate_particles ()
 
void initialize_particles ()
 
void advance_timestep ()
 
types::particle_index n_global_particles () const
 
void connect_to_signals (aspect::SimulatorSignals< dim > &signals)
 
unsigned int cell_weight (const typename parallel::distributed::Triangulation< dim >::cell_iterator &cell, const typename parallel::distributed::Triangulation< dim >::CellStatus status)
 
void update_particles ()
 
template<class Archive >
void serialize (Archive &ar, const unsigned int version)
 
virtual void save (std::ostringstream &os) const
 
virtual void load (std::istringstream &is)
 
virtual void parse_parameters (ParameterHandler &prm)
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()=default
 
virtual void initialize_simulator (const Simulator< dim > &simulator_object)
 
const Introspection< dim > & introspection () const
 
const Simulator< dim > & get_simulator () const
 
const Parameters< dim > & get_parameters () const
 
SimulatorSignals< dim > & get_signals () const
 
MPI_Comm get_mpi_communicator () const
 
TimerOutput & get_computing_timer () const
 
const ConditionalOStream & get_pcout () const
 
double get_time () const
 
double get_timestep () const
 
double get_old_timestep () const
 
unsigned int get_timestep_number () const
 
const TimeStepping::Manager< dim > & get_timestepping_manager () const
 
unsigned int get_nonlinear_iteration () const
 
const parallel::distributed::Triangulation< dim > & get_triangulation () const
 
double get_volume () const
 
const Mapping< dim > & get_mapping () const
 
std::string get_output_directory () const
 
bool include_adiabatic_heating () const
 
bool include_latent_heat () const
 
bool include_melt_transport () const
 
int get_stokes_velocity_degree () const
 
double get_adiabatic_surface_temperature () const
 
double get_surface_pressure () const
 
bool convert_output_to_years () const
 
unsigned int get_pre_refinement_step () const
 
unsigned int n_compositional_fields () const
 
double get_end_time () const
 
void get_refinement_criteria (Vector< float > &estimated_error_per_cell) const
 
void get_artificial_viscosity (Vector< float > &viscosity_per_cell, const bool skip_interior_cells=false) const
 
void get_artificial_viscosity_composition (Vector< float > &viscosity_per_cell, const unsigned int compositional_variable) const
 
const LinearAlgebra::BlockVectorget_current_linearization_point () const
 
const LinearAlgebra::BlockVectorget_solution () const
 
const LinearAlgebra::BlockVectorget_old_solution () const
 
const LinearAlgebra::BlockVectorget_old_old_solution () const
 
const LinearAlgebra::BlockVectorget_reaction_vector () const
 
const LinearAlgebra::BlockVectorget_mesh_velocity () const
 
const DoFHandler< dim > & get_dof_handler () const
 
const FiniteElement< dim > & get_fe () const
 
const LinearAlgebra::BlockSparseMatrixget_system_matrix () const
 
const LinearAlgebra::BlockSparseMatrixget_system_preconditioner_matrix () const
 
const MaterialModel::Interface< dim > & get_material_model () const
 
const GravityModel::Interface< dim > & get_gravity_model () const
 
const InitialTopographyModel::Interface< dim > & get_initial_topography_model () const
 
const GeometryModel::Interface< dim > & get_geometry_model () const
 
const AdiabaticConditions::Interface< dim > & get_adiabatic_conditions () const
 
bool has_boundary_temperature () const
 
const BoundaryTemperature::Manager< dim > & get_boundary_temperature_manager () const
 
const BoundaryHeatFlux::Interface< dim > & get_boundary_heat_flux () const
 
bool has_boundary_composition () const
 
const BoundaryComposition::Manager< dim > & get_boundary_composition_manager () const
 
const BoundaryTraction::Manager< dim > & get_boundary_traction_manager () const
 
std::shared_ptr< const InitialTemperature::Manager< dim > > get_initial_temperature_manager_pointer () const
 
const InitialTemperature::Manager< dim > & get_initial_temperature_manager () const
 
std::shared_ptr< const InitialComposition::Manager< dim > > get_initial_composition_manager_pointer () const
 
const InitialComposition::Manager< dim > & get_initial_composition_manager () const
 
const std::set< types::boundary_id > & get_fixed_temperature_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_heat_flux_boundary_indicators () const
 
const std::set< types::boundary_id > & get_fixed_composition_boundary_indicators () const
 
const std::set< types::boundary_id > & get_mesh_deformation_boundary_indicators () const
 
const BoundaryVelocity::Manager< dim > & get_boundary_velocity_manager () const
 
const HeatingModel::Manager< dim > & get_heating_model_manager () const
 
const MeshRefinement::Manager< dim > & get_mesh_refinement_manager () const
 
const MeltHandler< dim > & get_melt_handler () const
 
const VolumeOfFluidHandler< dim > & get_volume_of_fluid_handler () const
 
const NewtonHandler< dim > & get_newton_handler () const
 
const MeshDeformation::MeshDeformationHandler< dim > & get_mesh_deformation_handler () const
 
const LateralAveraging< dim > & get_lateral_averaging () const
 
const AffineConstraints< double > & get_current_constraints () const
 
bool simulator_is_past_initialization () const
 
double get_pressure_scaling () const
 
bool pressure_rhs_needs_compatibility_modification () const
 
bool model_has_prescribed_stokes_solution () const
 
TableHandler & get_statistics_object () const
 
const Postprocess::Manager< dim > & get_postprocess_manager () const
 
const Particle::World< dim > & get_particle_world () const
 
Particle::World< dim > & get_particle_world ()
 
bool is_stokes_matrix_free ()
 
const StokesMatrixFreeHandler< dim > & get_stokes_matrix_free () const
 
RotationProperties< dim > compute_net_angular_momentum (const bool use_constant_density, const LinearAlgebra::BlockVector &solution, const bool limit_to_top_faces=false) const
 

Static Public Member Functions

static void declare_parameters (ParameterHandler &prm)
 
- Static Public Member Functions inherited from aspect::SimulatorAccess< dim >
static void get_composition_values_at_q_point (const std::vector< std::vector< double >> &composition_values, const unsigned int q, std::vector< double > &composition_values_at_q_point)
 

Private Member Functions

std::map< types::subdomain_id, unsigned int > get_subdomain_id_to_neighbor_map () const
 
void apply_particle_per_cell_bounds ()
 
void advect_particles ()
 
void local_initialize_particles (const typename ParticleHandler< dim >::particle_iterator &begin_particle, const typename ParticleHandler< dim >::particle_iterator &end_particle)
 
void local_update_particles (const typename DoFHandler< dim >::active_cell_iterator &cell, const typename ParticleHandler< dim >::particle_iterator &begin_particle, const typename ParticleHandler< dim >::particle_iterator &end_particle)
 
void local_update_particles (const typename DoFHandler< dim >::active_cell_iterator &cell, const typename ParticleHandler< dim >::particle_iterator &begin_particle, const typename ParticleHandler< dim >::particle_iterator &end_particle, internal::SolutionEvaluators< dim > &evaluators)
 
void local_advect_particles (const typename DoFHandler< dim >::active_cell_iterator &cell, const typename ParticleHandler< dim >::particle_iterator &begin_particle, const typename ParticleHandler< dim >::particle_iterator &end_particle)
 
void local_advect_particles (const typename DoFHandler< dim >::active_cell_iterator &cell, const typename ParticleHandler< dim >::particle_iterator &begin_particle, const typename ParticleHandler< dim >::particle_iterator &end_particle, internal::SolutionEvaluators< dim > &evaluators)
 
void connect_particle_handler_signals (aspect::SimulatorSignals< dim > &signals, ParticleHandler< dim > &particle_handler, const bool connect_to_checkpoint_signals=true) const
 

Private Attributes

std::unique_ptr< Generator::Interface< dim > > generator
 
std::unique_ptr< Integrator::Interface< dim > > integrator
 
std::unique_ptr< Interpolator::Interface< dim > > interpolator
 
std::unique_ptr< Particles::ParticleHandler< dim > > particle_handler
 
Particles::ParticleHandler< dim > particle_handler_backup
 
std::unique_ptr< Property::Manager< dim > > property_manager
 
ParticleLoadBalancing::Kind particle_load_balancing
 
unsigned int min_particles_per_cell
 
unsigned int max_particles_per_cell
 
unsigned int particle_weight
 
bool update_ghost_particles
 

Detailed Description

template<int dim>
class aspect::Particle::World< dim >

This class manages the storage and handling of particles. It provides interfaces to generate and store particles, functions to initialize, update and advect them, and ways to retrieve information about the particles. The implementation of most of these methods is outsourced to different plugin systems, this class is mostly concerned with managing the interactions of the different systems with the code outside the particle world.

Definition at line 143 of file simulator_access.h.

Constructor & Destructor Documentation

§ World()

template<int dim>
aspect::Particle::World< dim >::World ( )

Default World constructor.

§ ~World()

template<int dim>
aspect::Particle::World< dim >::~World ( )
override

Default World destructor.

Member Function Documentation

§ initialize()

template<int dim>
void aspect::Particle::World< dim >::initialize ( )

Initialize the particle world.

§ get_property_manager()

template<int dim>
const Property::Manager<dim>& aspect::Particle::World< dim >::get_property_manager ( ) const

Get the particle property manager for this particle world.

Returns
The property manager for this world.

§ get_particle_handler() [1/2]

template<int dim>
const Particles::ParticleHandler<dim>& aspect::Particle::World< dim >::get_particle_handler ( ) const

Get the const particle handler for this particle world.

Returns
The particle handler for this world.

§ get_particle_handler() [2/2]

template<int dim>
Particles::ParticleHandler<dim>& aspect::Particle::World< dim >::get_particle_handler ( )

Get the particle handler for this particle world. There is no get_particles() function in the deal.II ParticleHandler, so we get and set the positions of the particles. These getter/setter functions are not const, and neither are the calling functions, but the existing get_particle_handler is. Therefore this non-const function is added.

Returns
The particle handler for this world.

§ copy_particle_handler()

template<int dim>
void aspect::Particle::World< dim >::copy_particle_handler ( const Particles::ParticleHandler< dim > &  from_particle_handler,
Particles::ParticleHandler< dim > &  to_particle_handler 
) const

Copy the state of particle handler from_particle_handler into the particle handler to_particle_handler. This will copy all particles and properties and leave to_particle_handler as an identical copy of from_particle_handler, assuming it was set up by this particle world class. This means we assume from_particle_handler uses the same triangulation and particle properties as are used in this model. Existing particles in to_particle_handler are deleted.

This function is expensive as it has to duplicate all data in from_particle_handler, and insert it into to_particle_handler, 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.

§ backup_particles()

template<int dim>
void aspect::Particle::World< dim >::backup_particles ( )

Stores a copy of the particle handler in particle_handler_backup. This copy can be used to restore the position and properties of the particles for example after advection solver iterations of the iterated advection scheme or when a timestep has to be repeated.

§ restore_particles()

template<int dim>
void aspect::Particle::World< dim >::restore_particles ( )

Restores the particle handler particle_handler based on the copy in particle_handler_backup. This restores the position and properties of the particles to those in the copy and can be used for example after advection solver iterations of the iterated advection scheme or when a timestep has to be repeated.

§ setup_initial_state()

template<int dim>
void aspect::Particle::World< dim >::setup_initial_state ( )

Do initial logic for handling pre-refinement steps

§ get_interpolator()

template<int dim>
const Interpolator::Interface<dim>& aspect::Particle::World< dim >::get_interpolator ( ) const

Get the particle interpolator for this particle world.

Returns
The interpolator for this world.

§ generate_particles()

template<int dim>
void aspect::Particle::World< dim >::generate_particles ( )

Initialize the particle properties.

§ initialize_particles()

template<int dim>
void aspect::Particle::World< dim >::initialize_particles ( )

Initialize the particle properties.

§ advance_timestep()

template<int dim>
void aspect::Particle::World< dim >::advance_timestep ( )

Advance particles by the old timestep using the current integration scheme. This accounts for the fact that the particles are actually still at their old positions and the current timestep length is already updated for the next step at the time this function is called.

§ n_global_particles()

template<int dim>
types::particle_index aspect::Particle::World< dim >::n_global_particles ( ) const

Return the total number of particles in the simulation. This function is useful for monitoring how many particles have been lost by falling out of the domain. Note that this function does not compute the number of particles, because that is an expensive global MPI operation. Instead it returns the number, which is updated internally every time it might change by a call to update_n_global_particles().

Returns
Total number of particles in simulation.

§ connect_to_signals()

template<int dim>
void aspect::Particle::World< dim >::connect_to_signals ( aspect::SimulatorSignals< dim > &  signals)

This callback function is registered within Simulator by the constructor of this class and will be called from Simulator during construction. It allows to attach slot functions to the provided SimulatorSignals. This world will register the register_store_callback_function() and register_load_callback_function() to the pre_refinement_store_user_data signal and the post_refinement_load_user_data signal respectively.

§ cell_weight()

template<int dim>
unsigned int aspect::Particle::World< dim >::cell_weight ( const typename parallel::distributed::Triangulation< dim >::cell_iterator &  cell,
const typename parallel::distributed::Triangulation< dim >::CellStatus  status 
)

Called by listener functions from Triangulation for every cell before a refinement step. A weight is attached to every cell depending on the number of contained particles.

§ update_particles()

template<int dim>
void aspect::Particle::World< dim >::update_particles ( )

Update the particle properties if necessary.

§ serialize()

template<int dim>
template<class Archive >
void aspect::Particle::World< dim >::serialize ( Archive &  ar,
const unsigned int  version 
)

Serialize the contents of this class.

Definition at line 545 of file world.h.

§ save()

template<int dim>
virtual void aspect::Particle::World< dim >::save ( std::ostringstream &  os) const
virtual

Save the state of the object.

§ load()

template<int dim>
virtual void aspect::Particle::World< dim >::load ( std::istringstream &  is)
virtual

Restore the state of the object.

§ declare_parameters()

template<int dim>
static void aspect::Particle::World< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this class takes through input files.

§ parse_parameters()

template<int dim>
virtual void aspect::Particle::World< dim >::parse_parameters ( ParameterHandler &  prm)
virtual

Read the parameters this class declares from the parameter file.

§ get_subdomain_id_to_neighbor_map()

template<int dim>
std::map<types::subdomain_id, unsigned int> aspect::Particle::World< dim >::get_subdomain_id_to_neighbor_map ( ) const
private

Get a map between subdomain id and the neighbor index. In other words the returned map answers the question: Given a subdomain id, which neighbor of the current processor's domain (in terms of a contiguous number from 0 to n_neighbors) owns this subdomain?

§ apply_particle_per_cell_bounds()

template<int dim>
void aspect::Particle::World< dim >::apply_particle_per_cell_bounds ( )
private

Apply the bounds for the maximum and minimum number of particles per cell, if the appropriate particle_load_balancing strategy has been selected.

§ advect_particles()

template<int dim>
void aspect::Particle::World< dim >::advect_particles ( )
private

Advect the particle positions by one integration step. Needs to be called until integrator->continue() returns false.

§ local_initialize_particles()

template<int dim>
void aspect::Particle::World< dim >::local_initialize_particles ( const typename ParticleHandler< dim >::particle_iterator &  begin_particle,
const typename ParticleHandler< dim >::particle_iterator &  end_particle 
)
private

Initialize the particle properties of one cell.

§ local_update_particles() [1/2]

template<int dim>
void aspect::Particle::World< dim >::local_update_particles ( const typename DoFHandler< dim >::active_cell_iterator &  cell,
const typename ParticleHandler< dim >::particle_iterator &  begin_particle,
const typename ParticleHandler< dim >::particle_iterator &  end_particle 
)
private

Update the particle properties of one cell.

§ local_update_particles() [2/2]

template<int dim>
void aspect::Particle::World< dim >::local_update_particles ( const typename DoFHandler< dim >::active_cell_iterator &  cell,
const typename ParticleHandler< dim >::particle_iterator &  begin_particle,
const typename ParticleHandler< dim >::particle_iterator &  end_particle,
internal::SolutionEvaluators< dim > &  evaluators 
)
private

Update the particle properties of one cell.

This version of the function above uses the deal.II class FEPointEvaluation, which allows for a faster evaluation of the solution under certain conditions. Check the place where this function is called for a list of conditions (they will change while FEPointEvaluation is generalized).

§ local_advect_particles() [1/2]

template<int dim>
void aspect::Particle::World< dim >::local_advect_particles ( const typename DoFHandler< dim >::active_cell_iterator &  cell,
const typename ParticleHandler< dim >::particle_iterator &  begin_particle,
const typename ParticleHandler< dim >::particle_iterator &  end_particle 
)
private

Advect the particles of one cell. Performs only one step for multi-step integrators. Needs to be called until integrator->continue() evaluates to false. Particles that moved out of their old cell during this advection step are removed from the local multimap and stored in particles_out_of_cell for further treatment (sorting them into the new cell).

§ local_advect_particles() [2/2]

template<int dim>
void aspect::Particle::World< dim >::local_advect_particles ( const typename DoFHandler< dim >::active_cell_iterator &  cell,
const typename ParticleHandler< dim >::particle_iterator &  begin_particle,
const typename ParticleHandler< dim >::particle_iterator &  end_particle,
internal::SolutionEvaluators< dim > &  evaluators 
)
private

Advect the particles of one cell. Performs only one step for multi-step integrators. Needs to be called until integrator->continue() evaluates to false. Particles that moved out of their old cell during this advection step are removed from the local multimap and stored in particles_out_of_cell for further treatment (sorting them into the new cell).

This version of the above function uses the deal.II class FEPointEvaluation, which allows for a faster evaluation of the solution under certain conditions. Check the place where this function is called for a list of conditions (they will change while FEPointEvaluation is generalized).

§ connect_particle_handler_signals()

template<int dim>
void aspect::Particle::World< dim >::connect_particle_handler_signals ( aspect::SimulatorSignals< dim > &  signals,
ParticleHandler< dim > &  particle_handler,
const bool  connect_to_checkpoint_signals = true 
) const
private

This function registers the necessary functions to the signals that the particle_handler needs to know about.

Member Data Documentation

§ generator

template<int dim>
std::unique_ptr<Generator::Interface<dim> > aspect::Particle::World< dim >::generator
private

Generation scheme for creating particles in this world

Definition at line 363 of file world.h.

§ integrator

template<int dim>
std::unique_ptr<Integrator::Interface<dim> > aspect::Particle::World< dim >::integrator
private

Integration scheme for moving particles in this world

Definition at line 368 of file world.h.

§ interpolator

template<int dim>
std::unique_ptr<Interpolator::Interface<dim> > aspect::Particle::World< dim >::interpolator
private

Interpolation scheme for moving particles in this world

Definition at line 373 of file world.h.

§ particle_handler

template<int dim>
std::unique_ptr<Particles::ParticleHandler<dim> > aspect::Particle::World< dim >::particle_handler
private

Particle handler object that is responsible for storing and managing the internal particle structures.

Definition at line 379 of file world.h.

§ particle_handler_backup

template<int dim>
Particles::ParticleHandler<dim> aspect::Particle::World< dim >::particle_handler_backup
private

Particle handler object that is responsible for storing and managing the internal particle structures before starting nonlinear iterations such that the particle position and properties can be restored after each outer advection iteration.

Definition at line 387 of file world.h.

§ property_manager

template<int dim>
std::unique_ptr<Property::Manager<dim> > aspect::Particle::World< dim >::property_manager
private

The property manager stores information about the additional particle properties and handles the initialization and update of these properties.

Definition at line 394 of file world.h.

§ particle_load_balancing

template<int dim>
ParticleLoadBalancing::Kind aspect::Particle::World< dim >::particle_load_balancing
private

Strategy for particle load balancing.

Definition at line 399 of file world.h.

§ min_particles_per_cell

template<int dim>
unsigned int aspect::Particle::World< dim >::min_particles_per_cell
private

Lower limit for particle number per cell. This limit is useful for adaptive meshes to prevent fine cells from being empty of particles. It will be checked and enforced after mesh refinement and after particle movement. If there are n_number_of_particles < min_particles_per_cell particles in one cell then min_particles_per_cell - n_number_of_particles particles are generated and randomly placed in this cell. If the particles carry properties the individual property plugins control how the properties of the new particles are initialized.

Definition at line 413 of file world.h.

§ max_particles_per_cell

template<int dim>
unsigned int aspect::Particle::World< dim >::max_particles_per_cell
private

Upper limit for particle number per cell. This limit is useful for adaptive meshes to prevent coarse cells from slowing down the whole model. It will be checked and enforced after mesh refinement, after MPI transfer of particles and after particle movement. If there are n_number_of_particles > max_particles_per_cell particles in one cell then n_number_of_particles - max_particles_per_cell particles in this cell are randomly chosen and destroyed.

Definition at line 426 of file world.h.

§ particle_weight

template<int dim>
unsigned int aspect::Particle::World< dim >::particle_weight
private

The computational cost of a single particle. This is an input parameter that is set during initialization and is only used if the particle load balancing strategy 'repartition' is used. This value determines how costly the computation of a single particle is compared to the computation of a whole cell, which is arbitrarily defined to represent a cost of 1000.

Definition at line 436 of file world.h.

§ update_ghost_particles

template<int dim>
bool aspect::Particle::World< dim >::update_ghost_particles
private

Some particle interpolation algorithms require knowledge about particles in neighboring cells. To allow this, particles in ghost cells need to be exchanged between the processes neighboring this cell. This parameter determines whether this transport is happening.

Definition at line 445 of file world.h.


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