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

Public Member Functions

 Manager ()
 
 ~Manager () override
 
void initialize ()
 
void initialize_one_particle (typename ParticleHandler< dim >::particle_iterator &particle) const
 
std::vector< double > initialize_late_particle (const Point< dim > &particle_location, const ParticleHandler< dim > &particle_handler, const Interpolator::Interface< dim > &interpolator, const typename parallel::distributed::Triangulation< dim >::active_cell_iterator &cell=typename parallel::distributed::Triangulation< dim >::active_cell_iterator()) const
 
void update_one_particle (typename ParticleHandler< dim >::particle_iterator &particle, const Vector< double > &solution, const std::vector< Tensor< 1, dim >> &gradients) const
 
UpdateTimeFlags need_update () const
 
UpdateFlags get_needed_update_flags () const
 
bool plugin_name_exists (const std::string &name) const
 
bool check_plugin_order (const std::string &first, const std::string &second) const
 
unsigned int get_plugin_index_by_name (const std::string &name) const
 
template<typename ParticlePropertyType , typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,ParticlePropertyType>::value>>
bool has_matching_property () const
 
template<typename ParticlePropertyType , typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,ParticlePropertyType>::value>>
const ParticlePropertyType & get_matching_property () const
 
unsigned int get_n_property_components () const
 
std::size_t get_particle_size () const
 
const ParticlePropertyInformationget_data_info () const
 
DEAL_II_DEPRECATED unsigned int get_property_component_by_name (const std::string &name) const
 
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 register_particle_property (const std::string &name, const std::string &description, void(*declare_parameters_function)(ParameterHandler &), std::unique_ptr< Property::Interface< dim >>(*factory_function)())
 
static void declare_parameters (ParameterHandler &prm)
 
static void write_plugin_graph (std::ostream &output_stream)
 
- 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 Attributes

std::vector< std::string > plugin_names
 
std::list< std::unique_ptr< Interface< dim > > > property_list
 
ParticlePropertyInformation property_information
 

Detailed Description

template<int dim>
class aspect::Particle::Property::Manager< dim >

Manager class of properties - This class sets the data of the collection of particles and updates it over time if requested by the user selected properties.

Definition at line 67 of file world.h.

Constructor & Destructor Documentation

§ Manager()

template<int dim>
aspect::Particle::Property::Manager< dim >::Manager ( )

Empty constructor for Manager

§ ~Manager()

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

Destructor for Manager

Member Function Documentation

§ initialize()

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

Initialization function. This function is called once at the beginning of the program after parse_parameters is run.

§ initialize_one_particle()

template<int dim>
void aspect::Particle::Property::Manager< dim >::initialize_one_particle ( typename ParticleHandler< dim >::particle_iterator &  particle) const

Initialization function for particle properties. This function is called once for each of the particles of a particle collection after it was created.

§ initialize_late_particle()

template<int dim>
std::vector<double> aspect::Particle::Property::Manager< dim >::initialize_late_particle ( const Point< dim > &  particle_location,
const ParticleHandler< dim > &  particle_handler,
const Interpolator::Interface< dim > &  interpolator,
const typename parallel::distributed::Triangulation< dim >::active_cell_iterator &  cell = typename parallel::distributed::Triangulation< dim >::active_cell_iterator() 
) const

Initialization function for particle properties. This function is called once for each of the particles of a particle collection that were created later than the initial particle generation.

§ update_one_particle()

template<int dim>
void aspect::Particle::Property::Manager< dim >::update_one_particle ( typename ParticleHandler< dim >::particle_iterator &  particle,
const Vector< double > &  solution,
const std::vector< Tensor< 1, dim >> &  gradients 
) const

Update function for particle properties. This function is called once every time step for every particle.

§ need_update()

template<int dim>
UpdateTimeFlags aspect::Particle::Property::Manager< dim >::need_update ( ) const

Returns an enum, which denotes at what time this class needs to update particle properties. The result of this class is a combination of the need_update() functions of all individual properties that are selected. More precise, it will choose to update the particle properties as often as the plugin that needs the most frequent update option requires. This saves considerable computation time, e.g. in cases when no plugin needs to update particle properties over time, because the solution does not need to be evaluated in this case.

§ get_needed_update_flags()

template<int dim>
UpdateFlags aspect::Particle::Property::Manager< dim >::get_needed_update_flags ( ) const

Return which data has to be provided to update all properties. Note that particle properties can only ask for update_default (no data), update_values (solution values), and update_gradients (solution gradients). All other update flags will have no effect.

§ plugin_name_exists()

template<int dim>
bool aspect::Particle::Property::Manager< dim >::plugin_name_exists ( const std::string &  name) const

Checks if the particle plugin specified by name exists in this model.

§ check_plugin_order()

template<int dim>
bool aspect::Particle::Property::Manager< dim >::check_plugin_order ( const std::string &  first,
const std::string &  second 
) const

Checks if the particle property plugin specified by first is executed before another particle property plugin specified by second.

Throws an assert when one of the plugin names does not exist. You can use the function plugin_name_exists() to check in advance whether a plugin exists

§ get_plugin_index_by_name()

template<int dim>
unsigned int aspect::Particle::Property::Manager< dim >::get_plugin_index_by_name ( const std::string &  name) const

Get the plugin index of the particle plugin specified by name.

§ has_matching_property()

template<int dim>
template<typename ParticlePropertyType , typename >
bool aspect::Particle::Property::Manager< dim >::has_matching_property ( ) const
inline

Go through the list of all particle properties that have been selected in the input file (and are consequently currently active) and return true if one of them has the desired type specified by the template argument.

This function can only be called if the given template type (the first template argument) is a class derived from the Interface class in this namespace.

Definition at line 805 of file interface.h.

§ get_matching_property()

template<int dim>
template<typename ParticlePropertyType , typename >
const ParticlePropertyType & aspect::Particle::Property::Manager< dim >::get_matching_property ( ) const
inline

Go through the list of all particle properties that have been selected in the input file (and are consequently currently active) and see if one of them has the type specified by the template argument or can be casted to that type. If so, return a reference to it. If no property is active that matches the given type, throw an exception.

This function can only be called if the given template type (the first template argument) is a class derived from the Interface class in this namespace.

Definition at line 819 of file interface.h.

§ get_n_property_components()

template<int dim>
unsigned int aspect::Particle::Property::Manager< dim >::get_n_property_components ( ) const

Get the number of components required to represent this particle's properties.

Returns
Number of doubles required to represent this particle's additional properties.

§ get_particle_size()

template<int dim>
std::size_t aspect::Particle::Property::Manager< dim >::get_particle_size ( ) const

Get the size in number of bytes required to represent this particle's properties for communication. This is essentially the space needed for the property components plus the space for the particle position and the space needed for its ID.

Returns
Number of bytes required to represent this particle.

§ get_data_info()

template<int dim>
const ParticlePropertyInformation& aspect::Particle::Property::Manager< dim >::get_data_info ( ) const

Get the names and number of components of particle properties.

Returns
A vector of pairs for each property name and the corresponding number of components attached to particles.

§ get_property_component_by_name()

template<int dim>
DEAL_II_DEPRECATED unsigned int aspect::Particle::Property::Manager< dim >::get_property_component_by_name ( const std::string &  name) const

Get the position of the property specified by name in the property vector of the particles.

Deprecated:
This function will be replaced by ParticlePropertyInformation::get_position_by_field_name(name)

§ register_particle_property()

template<int dim>
static void aspect::Particle::Property::Manager< dim >::register_particle_property ( const std::string &  name,
const std::string &  description,
void(*)(ParameterHandler &)  declare_parameters_function,
std::unique_ptr< Property::Interface< dim >>(*)()  factory_function 
)
static

A function that is used to register particle property objects in such a way that the Manager can deal with all of them without having to know them by name. This allows the files in which individual properties are implemented to register these properties, rather than also having to modify the Manager class by adding the new properties class.

Parameters
nameThe name under which this particle property is to be called in parameter files.
descriptionA text description of what this particle property does and that will be listed in the documentation of the parameter file.
declare_parameters_functionA pointer to a function that declares the parameters for this property.
factory_functionA pointer to a function that creates such a property object and returns a pointer to it.

§ declare_parameters()

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

Declare the parameters this class takes through input files.

§ parse_parameters()

template<int dim>
void aspect::Particle::Property::Manager< dim >::parse_parameters ( ParameterHandler &  prm)

Read the parameters this class declares from the parameter file.

§ write_plugin_graph()

template<int dim>
static void aspect::Particle::Property::Manager< dim >::write_plugin_graph ( std::ostream &  output_stream)
static

For the current plugin subsystem, write a connection graph of all of the plugins we know about, in the format that the programs dot and neato understand. This allows for a visualization of how all of the plugins that ASPECT knows about are interconnected, and connect to other parts of the ASPECT code.

Parameters
output_streamThe stream to write the output to.

Member Data Documentation

§ plugin_names

template<int dim>
std::vector<std::string> aspect::Particle::Property::Manager< dim >::plugin_names
private

Stores the names of the plugins which are present in the order they are executed.

Definition at line 783 of file interface.h.

§ property_list

template<int dim>
std::list<std::unique_ptr<Interface<dim> > > aspect::Particle::Property::Manager< dim >::property_list
private

A list of property objects that have been requested in the parameter file.

Definition at line 789 of file interface.h.

§ property_information

template<int dim>
ParticlePropertyInformation aspect::Particle::Property::Manager< dim >::property_information
private

A class that stores all information about the particle properties, their association with property plugins and their storage pattern.

Definition at line 795 of file interface.h.


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