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

Public Member Functions

 Diffusion ()
 
void initialize () override
 
void update () override
 
void compute_velocity_constraints_on_boundary (const DoFHandler< dim > &mesh_deformation_dof_handler, AffineConstraints< double > &mesh_velocity_constraints, const std::set< types::boundary_id > &boundary_id) const override
 
bool needs_surface_stabilization () const override
 
void parse_parameters (ParameterHandler &prm) override
 
- Public Member Functions inherited from aspect::MeshDeformation::Interface< dim >
virtual ~Interface ()=default
 
virtual Tensor< 1, dim > compute_initial_deformation_on_boundary (const types::boundary_id boundary_indicator, const Point< dim > &position) const
 
- 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::MeshDeformation::Interface< dim >
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

void diffuse_boundary (const DoFHandler< dim > &free_surface_dof_handler, const IndexSet &mesh_locally_owned, const IndexSet &mesh_locally_relevant, LinearAlgebra::Vector &output, const std::set< types::boundary_id > &boundary_id) const
 
void check_diffusion_time_step (const DoFHandler< dim > &mesh_deformation_dof_handler, const std::set< types::boundary_id > &boundary_ids) const
 

Private Attributes

double diffusivity
 
unsigned int timesteps_between_diffusion
 
bool apply_diffusion
 
std::set< types::boundary_id > additional_tangential_mesh_boundary_indicators
 

Detailed Description

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

A plugin that computes the deformation of surface vertices according to the solution of a dim-1 diffusion problem.

Definition at line 43 of file diffusion.h.

Constructor & Destructor Documentation

§ Diffusion()

template<int dim>
aspect::MeshDeformation::Diffusion< dim >::Diffusion ( )

Member Function Documentation

§ initialize()

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

Initialize function, which sets the start time and start timestep of diffusion.

Reimplemented from aspect::MeshDeformation::Interface< dim >.

§ update()

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

The update function sets the current time and determines whether diffusion should be applied in this timestep.

Reimplemented from aspect::MeshDeformation::Interface< dim >.

§ compute_velocity_constraints_on_boundary()

template<int dim>
void aspect::MeshDeformation::Diffusion< dim >::compute_velocity_constraints_on_boundary ( const DoFHandler< dim > &  mesh_deformation_dof_handler,
AffineConstraints< double > &  mesh_velocity_constraints,
const std::set< types::boundary_id > &  boundary_id 
) const
overridevirtual

A function that creates constraints for the velocity of certain mesh vertices (e.g. the surface vertices) for a specific boundary. The calling class will respect these constraints when computing the new vertex positions.

Reimplemented from aspect::MeshDeformation::Interface< dim >.

§ needs_surface_stabilization()

template<int dim>
bool aspect::MeshDeformation::Diffusion< dim >::needs_surface_stabilization ( ) const
overridevirtual

Returns whether or not the plugin requires surface stabilization

Reimplemented from aspect::MeshDeformation::Interface< dim >.

§ declare_parameters()

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

Declare parameters for the diffusion of the surface.

§ parse_parameters()

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

Parse parameters for the diffusion of the surface.

Reimplemented from aspect::MeshDeformation::Interface< dim >.

§ diffuse_boundary()

template<int dim>
void aspect::MeshDeformation::Diffusion< dim >::diffuse_boundary ( const DoFHandler< dim > &  free_surface_dof_handler,
const IndexSet &  mesh_locally_owned,
const IndexSet &  mesh_locally_relevant,
LinearAlgebra::Vector output,
const std::set< types::boundary_id > &  boundary_id 
) const
private

Compute the surface velocity from a difference in surface height given by the solution of the hillslope diffusion problem.

§ check_diffusion_time_step()

template<int dim>
void aspect::MeshDeformation::Diffusion< dim >::check_diffusion_time_step ( const DoFHandler< dim > &  mesh_deformation_dof_handler,
const std::set< types::boundary_id > &  boundary_ids 
) const
private

Check that the size of the next time step is not larger than the conduction timestep based on the mesh size and the diffusivity on each cell along the surface. The computed time step has to satisfy the CFL number chosen in the input parameter file on each cell of the mesh.

Member Data Documentation

§ diffusivity

template<int dim>
double aspect::MeshDeformation::Diffusion< dim >::diffusivity
private

The hillslope transport coefficient or diffusivity [m2/s] used in the hillslope diffusion of the deformed surface.

Definition at line 115 of file diffusion.h.

§ timesteps_between_diffusion

template<int dim>
unsigned int aspect::MeshDeformation::Diffusion< dim >::timesteps_between_diffusion
private

Maximum number of steps between the application of diffusion.

Definition at line 120 of file diffusion.h.

§ apply_diffusion

template<int dim>
bool aspect::MeshDeformation::Diffusion< dim >::apply_diffusion
private

Whether or not in the current timestep, diffusion should be applied.

Definition at line 126 of file diffusion.h.

§ additional_tangential_mesh_boundary_indicators

template<int dim>
std::set<types::boundary_id> aspect::MeshDeformation::Diffusion< dim >::additional_tangential_mesh_boundary_indicators
private

Boundaries along which the mesh is allowed to move tangentially despite of the Stokes velocity boundary conditions.

Definition at line 132 of file diffusion.h.


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