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

Public Member Functions

 MeshDeformationHandler (Simulator< dim > &simulator)
 
 ~MeshDeformationHandler () override
 
void initialize ()
 
void set_assemblers (const SimulatorAccess< dim > &simulator_access, aspect::Assemblers::Manager< dim > &assemblers) const
 
void update ()
 
void execute ()
 
void setup_dofs ()
 
void parse_parameters (ParameterHandler &prm)
 
const std::map< types::boundary_id, std::vector< std::string > > & get_active_mesh_deformation_names () const
 
const std::map< types::boundary_id, std::vector< std::unique_ptr< Interface< dim > > > > & get_active_mesh_deformation_models () const
 
const std::set< types::boundary_id > & get_active_mesh_deformation_boundary_indicators () const
 
const std::set< types::boundary_id > & get_boundary_indicators_requiring_stabilization () const
 
const std::set< types::boundary_id > & get_free_surface_boundary_indicators () const
 
double get_free_surface_theta () const
 
const LinearAlgebra::Vectorget_initial_topography () const
 
const LinearAlgebra::Vectorget_mesh_displacements () const
 
const DoFHandler< dim > & get_mesh_deformation_dof_handler () const
 
template<typename MeshDeformationType , typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,MeshDeformationType>::value>>
bool has_matching_mesh_deformation_object () const
 
template<typename MeshDeformationType , typename = typename std::enable_if_t<std::is_base_of<Interface<dim>,MeshDeformationType>::value>>
const MeshDeformationType & get_matching_mesh_deformation_object () const
 
const Mapping< dim > & get_level_mapping (const unsigned int level) const
 
 DeclException1 (ExcMeshDeformationNameNotFound, std::string,<< "Could not find entry <"<< arg1<< "> among the names of registered mesh deformation objects.")
 
- 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 void register_mesh_deformation (const std::string &name, const std::string &description, void(*declare_parameters_function)(ParameterHandler &), std::unique_ptr< Interface< dim >>(*factory_function)())
 
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 Member Functions

void make_initial_constraints ()
 
void make_constraints ()
 
void compute_mesh_displacements ()
 
void compute_mesh_displacements_gmg ()
 
void set_initial_topography ()
 
void interpolate_mesh_velocity ()
 
void update_multilevel_deformation ()
 

Private Attributes

Simulator< dim > & sim
 
const FESystem< dim > mesh_deformation_fe
 
DoFHandler< dim > mesh_deformation_dof_handler
 
LinearAlgebra::BlockVector mesh_velocity
 
LinearAlgebra::Vector mesh_displacements
 
LinearAlgebra::Vector old_mesh_displacements
 
LinearAlgebra::Vector initial_topography
 
LinearAlgebra::Vector fs_mesh_velocity
 
IndexSet mesh_locally_owned
 
IndexSet mesh_locally_relevant
 
AffineConstraints< double > mesh_velocity_constraints
 
AffineConstraints< double > mesh_vertex_constraints
 
std::map< types::boundary_id, std::vector< std::unique_ptr< Interface< dim > > > > mesh_deformation_objects
 
std::map< types::boundary_id, std::vector< std::string > > mesh_deformation_object_names
 
std::set< types::boundary_id > prescribed_mesh_deformation_boundary_indicators
 
std::set< types::boundary_id > tangential_mesh_deformation_boundary_indicators
 
std::set< types::boundary_id > zero_mesh_deformation_boundary_indicators
 
std::set< types::boundary_id > free_surface_boundary_indicators
 
std::set< types::boundary_id > boundary_indicators_requiring_stabilization
 
bool include_initial_topography
 
double surface_theta
 
MGLevelObject< std::unique_ptr< Mapping< dim > > > level_mappings
 
MGLevelObject<::LinearAlgebra::distributed::Vector< double > > level_displacements
 
MGTransferMF< dim, double > mg_transfer
 
MGConstrainedDoFs mg_constrained_dofs
 

Friends

class Simulator< dim >
 
class SimulatorAccess< dim >
 

Detailed Description

template<int dim>
class aspect::MeshDeformation::MeshDeformationHandler< dim >

The MeshDeformationHandler that handles the motion of the surface, the internal nodes and computes the Arbitrary-Lagrangian-Eulerian correction terms.

Definition at line 95 of file simulator.h.

Constructor & Destructor Documentation

§ MeshDeformationHandler()

template<int dim>
aspect::MeshDeformation::MeshDeformationHandler< dim >::MeshDeformationHandler ( Simulator< dim > &  simulator)

Initialize the mesh deformation handler, allowing it to read in relevant parameters as well as giving it a reference to the Simulator that owns it, since it needs to make fairly extensive changes to the internals of the simulator.

§ ~MeshDeformationHandler()

Destructor for the mesh deformation handler.

Member Function Documentation

§ initialize()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::initialize ( )

Initialization function of the MeshDeformationHandler.

The default implementation of this function does nothing.

§ set_assemblers()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::set_assemblers ( const SimulatorAccess< dim > &  simulator_access,
aspect::Assemblers::Manager< dim > &  assemblers 
) const

Called by Simulator::set_assemblers() to allow the FreeSurface plugin to register its assembler.

§ update()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::update ( )

Update function of the MeshDeformationHandler. This function allows the individual mesh deformation objects to update.

§ execute()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::execute ( )

The main execution step for the mesh deformation implementation. This computes the motion of the surface, moves the boundary nodes accordingly, redistributes the internal nodes in order to preserve mesh regularity, and calculates the Arbitrary-Lagrangian-Eulerian correction terms for advected quantities.

§ setup_dofs()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::setup_dofs ( )

Allocates and sets up the members of the MeshDeformationHandler. This is called by Simulator<dim>::setup_dofs()

§ declare_parameters()

template<int dim>
static void aspect::MeshDeformation::MeshDeformationHandler< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare parameters for the mesh deformation handling.

§ parse_parameters()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::parse_parameters ( ParameterHandler &  prm)

Parse parameters for the mesh deformation handling.

§ register_mesh_deformation()

template<int dim>
static void aspect::MeshDeformation::MeshDeformationHandler< dim >::register_mesh_deformation ( 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 mesh deformation 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 initial mesh deformation plugin class.

Parameters
nameA string that identifies the mesh deformation 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 mesh deformation model wants to read from input files.
factory_functionA pointer to a function that can create an object of this mesh deformation model.

§ get_active_mesh_deformation_names()

template<int dim>
const std::map<types::boundary_id, std::vector<std::string> >& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_active_mesh_deformation_names ( ) const

Return a map of boundary indicators to the names of all mesh deformation models currently used in the computation, as specified in the input file.

§ get_active_mesh_deformation_models()

template<int dim>
const std::map<types::boundary_id,std::vector<std::unique_ptr<Interface<dim> > > >& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_active_mesh_deformation_models ( ) const

Return a map of boundary indicators to vectors of pointers to all mesh deformation models currently used in the computation, as specified in the input file.

§ get_active_mesh_deformation_boundary_indicators()

template<int dim>
const std::set<types::boundary_id>& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_active_mesh_deformation_boundary_indicators ( ) const

Return a set of all the indicators of boundaries with mesh deformation objects on them.

§ get_boundary_indicators_requiring_stabilization()

template<int dim>
const std::set<types::boundary_id>& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_boundary_indicators_requiring_stabilization ( ) const

Return a set of all the indicators of boundaries that require surface stabilization.

§ get_free_surface_boundary_indicators()

template<int dim>
const std::set<types::boundary_id>& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_free_surface_boundary_indicators ( ) const

Return the boundary id of the surface that has a free surface mesh deformation object. If no free surface is used, an empty set is returned.

§ get_free_surface_theta()

template<int dim>
double aspect::MeshDeformation::MeshDeformationHandler< dim >::get_free_surface_theta ( ) const

Return the stabilization parameter for the free surface.

§ get_initial_topography()

template<int dim>
const LinearAlgebra::Vector& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_initial_topography ( ) const

Return the initial topography stored on the Q1 finite element that describes the mesh geometry. Note that a topography is set for all mesh nodes, but only the values of surface boundary nodes are correct. The internal nodes get the same initial topography as the corresponding surface node. In other words, there is no decrease of the initial topography with depth. However, only the topography stored at the surface nodes is taken into account in the diffusion plugin that uses this function. TODO Once all initial_topography is prescribed through initial_mesh_deformation, this function can be removed.

§ get_mesh_displacements()

template<int dim>
const LinearAlgebra::Vector& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_mesh_displacements ( ) const

Return the mesh displacements based on the mesh deformation DoFHandler you can access with get_mesh_deformation_dof_handler().

§ get_mesh_deformation_dof_handler()

template<int dim>
const DoFHandler<dim>& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_mesh_deformation_dof_handler ( ) const

Return the DoFHandler used to represent the mesh deformation space.

§ has_matching_mesh_deformation_object()

template<int dim>
template<typename MeshDeformationType , typename >
bool aspect::MeshDeformation::MeshDeformationHandler< dim >::has_matching_mesh_deformation_object ( ) const
inline

Go through the list of all mesh deformation objects 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 625 of file interface.h.

§ get_matching_mesh_deformation_object()

template<int dim>
template<typename MeshDeformationType , typename >
const MeshDeformationType & aspect::MeshDeformation::MeshDeformationHandler< dim >::get_matching_mesh_deformation_object ( ) const
inline

Go through the list of all mesh deformation objects 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 mesh deformation object 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 641 of file interface.h.

References aspect::MeshDeformation::get_valid_model_names_pattern().

§ get_level_mapping()

template<int dim>
const Mapping<dim>& aspect::MeshDeformation::MeshDeformationHandler< dim >::get_level_mapping ( const unsigned int  level) const

If multilevel solvers are used, we need a mapping on each multigrid level. These are automatically updated by this handler class and can be accessed with this method.

§ write_plugin_graph()

template<int dim>
static void aspect::MeshDeformation::MeshDeformationHandler< 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::MeshDeformation::MeshDeformationHandler< dim >::DeclException1 ( ExcMeshDeformationNameNotFound  ,
std::string   
)

Exception.

§ make_initial_constraints()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::make_initial_constraints ( )
private

Compute the initial constraints for the mesh displacement on the boundaries of the domain. This is used on the mesh deformation boundaries to describe a displacement (initial topography) to be used during the simulation. The displacement is given by the active deformation plugins.

§ make_constraints()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::make_constraints ( )
private

Compute the constraints for the mesh velocity on the boundaries of the domain. On the mesh deformation boundaries, the velocity is given by the active deformation plugins.

Velocities on free-slip boundaries are constrained to be tangential to those boundaries. Velocities on no-slip boundaries are set to be zero. If a no-slip boundary is marked as additional tangential, then velocities are constrained as tangential.

§ compute_mesh_displacements()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::compute_mesh_displacements ( )
private

Solve vector Laplacian equation for internal mesh displacements and update the current displacement vector based on the solution.

§ compute_mesh_displacements_gmg()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::compute_mesh_displacements_gmg ( )
private

Solve vector Laplacian equation using GMG for internal mesh displacements and update the current displacement vector based on the solution.

§ set_initial_topography()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::set_initial_topography ( )
private

Set up the vector with initial displacements of the mesh due to the initial topography, as supplied by the initial topography plugin based on the surface coordinates of the mesh nodes. We set all entries to the initial topography based on its surface coordinates, i.e. the initial topography is not corrected for depth from the surface as it is for the initial mesh deformation. TODO this is ok for now, because the surface diffusion plugin only cares about the initial topography at the surface, but it would be more correct if it sets the initial topography to the actual initial distortion of the mesh cells. When all initial_topography plugins are converted to the new initial_mesh_deformation functionality, this function can be removed.

§ interpolate_mesh_velocity()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::interpolate_mesh_velocity ( )
private

Calculate the velocity of the mesh for ALE corrections.

§ update_multilevel_deformation()

template<int dim>
void aspect::MeshDeformation::MeshDeformationHandler< dim >::update_multilevel_deformation ( )
private

Update the mesh deformation for the multigrid levels.

Friends And Related Function Documentation

§ Simulator< dim >

template<int dim>
friend class Simulator< dim >
friend

Definition at line 615 of file interface.h.

§ SimulatorAccess< dim >

template<int dim>
friend class SimulatorAccess< dim >
friend

Definition at line 616 of file interface.h.

Member Data Documentation

§ sim

template<int dim>
Simulator<dim>& aspect::MeshDeformation::MeshDeformationHandler< dim >::sim
private

Reference to the Simulator object to which a MeshDeformationHandler instance belongs.

Definition at line 465 of file interface.h.

§ mesh_deformation_fe

template<int dim>
const FESystem<dim> aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_deformation_fe
private

Finite element for the mesh deformation implementation, which is used for tracking mesh deformation.

Definition at line 471 of file interface.h.

§ mesh_deformation_dof_handler

template<int dim>
DoFHandler<dim> aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_deformation_dof_handler
private

DoFHandler for the mesh deformation implementation.

Definition at line 476 of file interface.h.

§ mesh_velocity

template<int dim>
LinearAlgebra::BlockVector aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_velocity
private

BlockVector which stores the mesh velocity in the Stokes finite element space. This is used for ALE corrections.

Definition at line 483 of file interface.h.

§ mesh_displacements

template<int dim>
LinearAlgebra::Vector aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_displacements
private

Vector for storing the positions of the mesh vertices. This is used for calculating the mapping from the reference cell to the position of the cell in the deformed mesh. This must be redistributed upon mesh refinement.

Definition at line 491 of file interface.h.

§ old_mesh_displacements

template<int dim>
LinearAlgebra::Vector aspect::MeshDeformation::MeshDeformationHandler< dim >::old_mesh_displacements
private

mesh_displacements from the last time step.

Definition at line 496 of file interface.h.

§ initial_topography

template<int dim>
LinearAlgebra::Vector aspect::MeshDeformation::MeshDeformationHandler< dim >::initial_topography
private

Vector for storing the positions of the mesh vertices at the initial timestep. This must be redistributed upon mesh refinement. We need to store the initial topography because it is not taken into account into the mesh displacements used by the MappingQ1Eulerian. The current mesh displacements plus the initial topography provide the actual topography at any time.

Definition at line 506 of file interface.h.

§ fs_mesh_velocity

template<int dim>
LinearAlgebra::Vector aspect::MeshDeformation::MeshDeformationHandler< dim >::fs_mesh_velocity
private

Vector for storing the mesh velocity in the mesh deformation finite element space, which is, in general, not the same finite element space as the Stokes system. This is used for interpolating the mesh velocity in the mesh deformation finite element space onto the velocity in the Stokes finite element space, which is then used for making the ALE correction in the advection equations.

Definition at line 516 of file interface.h.

§ mesh_locally_owned

template<int dim>
IndexSet aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_locally_owned
private

IndexSet for the locally owned DoFs for the mesh system

Definition at line 521 of file interface.h.

§ mesh_locally_relevant

template<int dim>
IndexSet aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_locally_relevant
private

IndexSet for the locally relevant DoFs for the mesh system

Definition at line 526 of file interface.h.

§ mesh_velocity_constraints

template<int dim>
AffineConstraints<double> aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_velocity_constraints
private

Storage for the mesh velocity constraints for solving the elliptic problem.

Definition at line 532 of file interface.h.

§ mesh_vertex_constraints

template<int dim>
AffineConstraints<double> aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_vertex_constraints
private

Storage for the mesh vertex constraints for keeping the mesh conforming upon redistribution.

Definition at line 538 of file interface.h.

§ mesh_deformation_objects

template<int dim>
std::map<types::boundary_id,std::vector<std::unique_ptr<Interface<dim> > > > aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_deformation_objects
private

A map of boundary ids to mesh deformation objects that have been requested in the parameter file.

Definition at line 544 of file interface.h.

§ mesh_deformation_object_names

template<int dim>
std::map<types::boundary_id, std::vector<std::string> > aspect::MeshDeformation::MeshDeformationHandler< dim >::mesh_deformation_object_names
private

Map from boundary id to a vector of names representing mesh deformation objects.

Definition at line 550 of file interface.h.

§ prescribed_mesh_deformation_boundary_indicators

template<int dim>
std::set<types::boundary_id> aspect::MeshDeformation::MeshDeformationHandler< dim >::prescribed_mesh_deformation_boundary_indicators
private

The set of boundary indicators for which mesh deformation objects are set and that therefore can deform over time as prescribed in the mesh_deformation_objects.

Definition at line 557 of file interface.h.

§ tangential_mesh_deformation_boundary_indicators

template<int dim>
std::set<types::boundary_id> aspect::MeshDeformation::MeshDeformationHandler< dim >::tangential_mesh_deformation_boundary_indicators
private

A set of boundary indicators that denote those boundaries that are allowed to move their mesh tangential to the boundary.

Definition at line 563 of file interface.h.

§ zero_mesh_deformation_boundary_indicators

template<int dim>
std::set<types::boundary_id> aspect::MeshDeformation::MeshDeformationHandler< dim >::zero_mesh_deformation_boundary_indicators
private

A set of boundary indicators, on which mesh deformation is prescribed to be zero (fixed boundaries that never move). All boundaries except those in prescribed_mesh_deformation_boundary_indicators and tangential_mesh_deformation_boundary_indicators are in this set.

Definition at line 571 of file interface.h.

§ free_surface_boundary_indicators

template<int dim>
std::set<types::boundary_id> aspect::MeshDeformation::MeshDeformationHandler< dim >::free_surface_boundary_indicators
private

The boundary indicator(s) of the free surface(s). This is the subset of prescribed_mesh_deformation_boundary_indicators for which the 'free surface' plugin was selected.

Definition at line 578 of file interface.h.

§ boundary_indicators_requiring_stabilization

template<int dim>
std::set<types::boundary_id> aspect::MeshDeformation::MeshDeformationHandler< dim >::boundary_indicators_requiring_stabilization
private

The set of boundary indicators for which the mesh deformation objects need surface stabilization.

Definition at line 584 of file interface.h.

§ include_initial_topography

template<int dim>
bool aspect::MeshDeformation::MeshDeformationHandler< dim >::include_initial_topography
private

Definition at line 586 of file interface.h.

§ surface_theta

template<int dim>
double aspect::MeshDeformation::MeshDeformationHandler< dim >::surface_theta
private

Stabilization parameter for the free surface. Should be between zero and one. A value of zero means no stabilization. See Kaus et. al. 2010 for more details.

Definition at line 593 of file interface.h.

§ level_mappings

template<int dim>
MGLevelObject<std::unique_ptr<Mapping<dim> > > aspect::MeshDeformation::MeshDeformationHandler< dim >::level_mappings
private

If required, store a mapping for each multigrid level.

Definition at line 598 of file interface.h.

§ level_displacements

template<int dim>
MGLevelObject<::LinearAlgebra::distributed::Vector<double> > aspect::MeshDeformation::MeshDeformationHandler< dim >::level_displacements
private

One vector on each multigrid level for the mesh displacement used in the mapping.

Definition at line 603 of file interface.h.

§ mg_transfer

template<int dim>
MGTransferMF<dim, double> aspect::MeshDeformation::MeshDeformationHandler< dim >::mg_transfer
private

Multigrid transfer operator for the displacements

Definition at line 608 of file interface.h.

§ mg_constrained_dofs

template<int dim>
MGConstrainedDoFs aspect::MeshDeformation::MeshDeformationHandler< dim >::mg_constrained_dofs
private

Multigrid level constraints for the displacements

Definition at line 613 of file interface.h.


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