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

Public Member Functions

Tensor< 1, dim > boundary_velocity (const types::boundary_id boundary_indicator, const Point< dim > &position) const
 
DEAL_II_DEPRECATED const std::map< types::boundary_id, std::pair< std::string, std::vector< std::string > > > & get_active_boundary_velocity_names () const
 
DEAL_II_DEPRECATED const std::map< types::boundary_id, std::vector< std::unique_ptr< BoundaryVelocity::Interface< dim > > > > & get_active_boundary_velocity_conditions () const
 
const std::vector< types::boundary_id > & get_active_plugin_boundary_indicators () const
 
ComponentMask get_component_mask (const types::boundary_id boundary_id) const
 
const std::set< types::boundary_id > & get_prescribed_boundary_velocity_indicators () const
 
const std::set< types::boundary_id > & get_zero_boundary_velocity_indicators () const
 
const std::set< types::boundary_id > & get_tangential_boundary_velocity_indicators () const
 
void parse_parameters (ParameterHandler &prm) override
 
template<typename BoundaryVelocityType , typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryVelocityType>::value>>
DEAL_II_DEPRECATED bool has_matching_boundary_velocity_model () const
 
template<typename BoundaryVelocityType , typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryVelocityType>::value>>
DEAL_II_DEPRECATED const BoundaryVelocityType & get_matching_boundary_velocity_model () const
 
 DeclException1 (ExcBoundaryVelocityNameNotFound, std::string,<< "Could not find entry <"<< arg1<< "> among the names of registered boundary velocity objects.")
 
template<typename BoundaryVelocityType , typename >
bool has_matching_boundary_velocity_model () const
 
template<typename BoundaryVelocityType , typename >
const BoundaryVelocityType & get_matching_boundary_velocity_model () const
 
- Public Member Functions inherited from aspect::Plugins::ManagerBase< Interface< dim > >
 ~ManagerBase () override
 
void update () override
 
bool has_matching_active_plugin () const
 
const PluginType & get_matching_active_plugin () const
 
const std::list< std::unique_ptr< Interface< dim > > > & get_active_plugins () const
 
const std::vector< std::string > & get_active_plugin_names () const
 
- Public Member Functions inherited from aspect::Plugins::InterfaceBase
virtual ~InterfaceBase ()=default
 
virtual void initialize ()
 
- 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
 
unsigned int n_particle_managers () const
 
const Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index) const
 
Particle::Manager< dim > & get_particle_manager (const unsigned int particle_manager_index)
 
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_boundary_velocity (const std::string &name, const std::string &description, void(*declare_parameters_function)(ParameterHandler &), std::unique_ptr< 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::Plugins::InterfaceBase
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 Attributes

std::vector< types::boundary_id > boundary_indicators
 
std::vector< ComponentMask > component_masks
 
std::map< types::boundary_id, std::vector< std::unique_ptr< BoundaryVelocity::Interface< dim > > > > boundary_velocity_objects
 
std::map< types::boundary_id, std::pair< std::string, std::vector< std::string > > > boundary_velocity_indicators
 
std::set< types::boundary_id > prescribed_velocity_boundary_indicators
 
std::set< types::boundary_id > zero_velocity_boundary_indicators
 
std::set< types::boundary_id > tangential_velocity_boundary_indicators
 

Additional Inherited Members

- Protected Attributes inherited from aspect::Plugins::ManagerBase< Interface< dim > >
std::list< std::unique_ptr< Interface< dim > > > plugin_objects
 
std::vector< std::string > plugin_names
 

Detailed Description

template<int dim>
class aspect::BoundaryVelocity::Manager< dim >

A class that manages all boundary velocity objects.

Definition at line 97 of file simulator_access.h.

Member Function Documentation

§ boundary_velocity()

template<int dim>
Tensor<1,dim> aspect::BoundaryVelocity::Manager< dim >::boundary_velocity ( const types::boundary_id  boundary_indicator,
const Point< dim > &  position 
) const

A function that calls the boundary_velocity functions of all the individual boundary velocity objects and uses the stored operators to combine them.

§ register_boundary_velocity()

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

A function that is used to register boundary velocity 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 plugins are implemented to register these plugins, rather than also having to modify the Manager class by adding the new boundary velocity plugin class.

Parameters
nameA string that identifies the boundary velocity model
descriptionA text description of what this model does and that will be listed in the documentation of the parameter file.
declare_parameters_functionA pointer to a function that can be used to declare the parameters that this boundary velocity model wants to read from input files.
factory_functionA pointer to a function that can create an object of this boundary velocity model.

§ get_active_boundary_velocity_names()

template<int dim>
DEAL_II_DEPRECATED const std::map<types::boundary_id, std::pair<std::string,std::vector<std::string> > >& aspect::BoundaryVelocity::Manager< dim >::get_active_boundary_velocity_names ( ) const

Return the names of all prescribed boundary velocity models currently used in the computation as specified in the input file. The function returns a map between a boundary identifier and a pair. The first part of the pair is a string that represents the prescribed velocity components on this boundary (e.g. y, xz, or xyz) and the second part is a vector of strings that represent the names of boundary velocity plugins for this boundary. If there are no prescribed boundary velocity plugins for a particular boundary, this boundary identifier will not appear in the map.

Deprecated:
This function will be removed. Use the function get_active_plugin_names() of the base class ManagerBase instead.

§ get_active_boundary_velocity_conditions()

template<int dim>
DEAL_II_DEPRECATED const std::map<types::boundary_id,std::vector<std::unique_ptr<BoundaryVelocity::Interface<dim> > > >& aspect::BoundaryVelocity::Manager< dim >::get_active_boundary_velocity_conditions ( ) const

Return pointers to all boundary velocity models currently used in the computation, as specified in the input file. The function returns a map between a boundary identifier and a vector of unique pointers that represent the names of prescribed velocity boundary models for this boundary. If there are no prescribed boundary velocity plugins for a particular boundary this boundary identifier will not appear in the map.

Deprecated:
This function has been removed. Use the function get_active_plugins() of the base class ManagerBase instead.

§ get_active_plugin_boundary_indicators()

template<int dim>
const std::vector<types::boundary_id>& aspect::BoundaryVelocity::Manager< dim >::get_active_plugin_boundary_indicators ( ) const

Return a list of boundary indicators that indicate for each active plugin which boundary id it is responsible for. The list of active plugins can be requested by calling get_active_plugins().

§ get_component_mask()

template<int dim>
ComponentMask aspect::BoundaryVelocity::Manager< dim >::get_component_mask ( const types::boundary_id  boundary_id) const

Return a component mask that indicates for the given boundary_id which velocity components are prescribed by this manager class. All plugins that are responsible for this boundary use the same component mask. The list of plugin objects can be requested by calling get_active_plugins() and the list of boundaries they are responsible for is returned by get_active_plugin_boundary_indicators().

§ get_prescribed_boundary_velocity_indicators()

template<int dim>
const std::set<types::boundary_id>& aspect::BoundaryVelocity::Manager< dim >::get_prescribed_boundary_velocity_indicators ( ) const

Return a set of boundary indicators for which boundary velocities are prescribed.

§ get_zero_boundary_velocity_indicators()

template<int dim>
const std::set<types::boundary_id>& aspect::BoundaryVelocity::Manager< dim >::get_zero_boundary_velocity_indicators ( ) const

Return a list of boundary ids on which the velocity is prescribed to be zero (no-slip).

§ get_tangential_boundary_velocity_indicators()

template<int dim>
const std::set<types::boundary_id>& aspect::BoundaryVelocity::Manager< dim >::get_tangential_boundary_velocity_indicators ( ) const

Return a list of boundary ids on which the velocity is prescribed to be tangential (free-slip).

§ declare_parameters()

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

Declare the parameters of all known boundary velocity plugins, as well as the ones this class has itself.

§ parse_parameters()

template<int dim>
void aspect::BoundaryVelocity::Manager< dim >::parse_parameters ( ParameterHandler &  prm)
overridevirtual

Read the parameters this class declares from the parameter file. This determines which boundary velocity objects will be created; then let these objects read their parameters as well.

Reimplemented from aspect::Plugins::InterfaceBase.

§ has_matching_boundary_velocity_model() [1/2]

template<int dim>
template<typename BoundaryVelocityType , typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryVelocityType>::value>>
DEAL_II_DEPRECATED bool aspect::BoundaryVelocity::Manager< dim >::has_matching_boundary_velocity_model ( ) const

Go through the list of all boundary velocity models 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.

Deprecated:
Instead of this function, use the Plugins::ManagerBase::has_matching_active_plugin() and Plugins::ManagerBase::get_matching_active_plugin() functions of the base class of the current class.

§ get_matching_boundary_velocity_model() [1/2]

template<int dim>
template<typename BoundaryVelocityType , typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,BoundaryVelocityType>::value>>
DEAL_II_DEPRECATED const BoundaryVelocityType& aspect::BoundaryVelocity::Manager< dim >::get_matching_boundary_velocity_model ( ) const

Go through the list of all boundary velocity models 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 cast to that type. If so, return a reference to it. If no boundary velocity model 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.

Deprecated:
Instead of this function, use the Plugins::ManagerBase::has_matching_active_plugin() and Plugins::ManagerBase::get_matching_active_plugin() functions of the base class of the current class.

§ write_plugin_graph()

template<int dim>
static void aspect::BoundaryVelocity::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.

§ DeclException1()

template<int dim>
aspect::BoundaryVelocity::Manager< dim >::DeclException1 ( ExcBoundaryVelocityNameNotFound  ,
std::string   
)

Exception.

§ has_matching_boundary_velocity_model() [2/2]

template<int dim>
template<typename BoundaryVelocityType , typename >
bool aspect::BoundaryVelocity::Manager< dim >::has_matching_boundary_velocity_model ( ) const
inline

Definition at line 344 of file interface.h.

§ get_matching_boundary_velocity_model() [2/2]

template<int dim>
template<typename BoundaryVelocityType , typename >
const BoundaryVelocityType& aspect::BoundaryVelocity::Manager< dim >::get_matching_boundary_velocity_model ( ) const
inline

Member Data Documentation

§ boundary_indicators

template<int dim>
std::vector<types::boundary_id> aspect::BoundaryVelocity::Manager< dim >::boundary_indicators
private

A list of boundary indicators that indicate for each plugin in the list of plugin_objects which boundary id it is responsible for. By default each plugin is active for all boundaries, but this list can be modified by derived classes to limit the application of plugins to specific boundaries.

Definition at line 284 of file interface.h.

§ component_masks

template<int dim>
std::vector<ComponentMask> aspect::BoundaryVelocity::Manager< dim >::component_masks
private

A list of boundary indicators that indicate for each plugin in the list of plugin_objects which components it is responsible for. By default each plugin is active for all components, but this list can be modified by derived classes to limit the application of plugins to specific boundaries.

Definition at line 294 of file interface.h.

§ boundary_velocity_objects

template<int dim>
std::map<types::boundary_id,std::vector<std::unique_ptr<BoundaryVelocity::Interface<dim> > > > aspect::BoundaryVelocity::Manager< dim >::boundary_velocity_objects
private

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

Deprecated:
This variable is no longer used, but needed to issue a proper error message in the function get_active_boundary_velocity_conditions().

Definition at line 303 of file interface.h.

§ boundary_velocity_indicators

template<int dim>
std::map<types::boundary_id, std::pair<std::string,std::vector<std::string> > > aspect::BoundaryVelocity::Manager< dim >::boundary_velocity_indicators
private

Map from boundary id to a pair ("components", list of "velocity boundary type"), where components is of the format "[x][y][z]" and the velocity type is mapped to one of the plugins of velocity boundary conditions (e.g. "function"). If the components string is empty, it is assumed the plugins are used for all components.

Deprecated:
Remove this variable when the deprecated functions get_active_boundary_velocity_names and get_active_boundary_velocity_conditions are removed. Use the base class variable plugin_names instead.

Definition at line 318 of file interface.h.

§ prescribed_velocity_boundary_indicators

template<int dim>
std::set<types::boundary_id> aspect::BoundaryVelocity::Manager< dim >::prescribed_velocity_boundary_indicators
private

A set of boundary indicators, on which velocities are prescribed.

Definition at line 323 of file interface.h.

§ zero_velocity_boundary_indicators

template<int dim>
std::set<types::boundary_id> aspect::BoundaryVelocity::Manager< dim >::zero_velocity_boundary_indicators
private

A set of boundary indicators, on which velocities are prescribed to zero (no-slip).

Definition at line 329 of file interface.h.

§ tangential_velocity_boundary_indicators

template<int dim>
std::set<types::boundary_id> aspect::BoundaryVelocity::Manager< dim >::tangential_velocity_boundary_indicators
private

A set of boundary indicators, on which velocities are prescribed to be tangential (free-slip).

Definition at line 335 of file interface.h.


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