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

Public Member Functions

void parse_parameters (ParameterHandler &prm)
 
void create_elastic_outputs (MaterialModel::MaterialModelOutputs< dim > &out) const
 
void fill_elastic_outputs (const MaterialModel::MaterialModelInputs< dim > &in, const std::vector< double > &average_elastic_shear_moduli, MaterialModel::MaterialModelOutputs< dim > &out) const
 
void fill_reaction_outputs (const MaterialModel::MaterialModelInputs< dim > &in, const std::vector< double > &average_elastic_shear_moduli, MaterialModel::MaterialModelOutputs< dim > &out) const
 
const std::vector< double > & get_elastic_shear_moduli () const
 
double calculate_elastic_viscosity (const double shear_modulus) const
 
double calculate_viscoelastic_viscosity (const double viscosity, const double shear_modulus) const
 
SymmetricTensor< 2, dim > calculate_viscoelastic_strain_rate (const SymmetricTensor< 2, dim > &strain_rate, const SymmetricTensor< 2, dim > &stored_stress, const double shear_modulus) const
 
double elastic_timestep () 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::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

double elastic_damper_viscosity
 
std::vector< double > elastic_shear_moduli
 
bool use_fixed_elastic_time_step
 
double fixed_elastic_time_step
 
double stabilization_time_scale_factor
 
std::unique_ptr< FEPointEvaluation< dim, dim > > evaluator
 

Detailed Description

template<int dim>
class aspect::MaterialModel::Rheology::Elasticity< dim >

Definition at line 62 of file elasticity.h.

Member Function Documentation

§ declare_parameters()

template<int dim>
static void aspect::MaterialModel::Rheology::Elasticity< dim >::declare_parameters ( ParameterHandler &  prm)
static

Declare the parameters this function takes through input files.

§ parse_parameters()

template<int dim>
void aspect::MaterialModel::Rheology::Elasticity< dim >::parse_parameters ( ParameterHandler &  prm)

Read the parameters from the parameter file.

§ create_elastic_outputs()

template<int dim>
void aspect::MaterialModel::Rheology::Elasticity< dim >::create_elastic_outputs ( MaterialModel::MaterialModelOutputs< dim > &  out) const

Create the additional material model outputs object that contains the elastic shear moduli.

§ fill_elastic_outputs()

template<int dim>
void aspect::MaterialModel::Rheology::Elasticity< dim >::fill_elastic_outputs ( const MaterialModel::MaterialModelInputs< dim > &  in,
const std::vector< double > &  average_elastic_shear_moduli,
MaterialModel::MaterialModelOutputs< dim > &  out 
) const

Given the stress of the previous time step in the material model inputs in, the elastic shear moduli average_elastic_shear_moduli a each point, and the (viscous) viscosities given in the material model outputs object out, fill an additional material model outputs objects with the elastic force terms.

§ fill_reaction_outputs()

template<int dim>
void aspect::MaterialModel::Rheology::Elasticity< dim >::fill_reaction_outputs ( const MaterialModel::MaterialModelInputs< dim > &  in,
const std::vector< double > &  average_elastic_shear_moduli,
MaterialModel::MaterialModelOutputs< dim > &  out 
) const

Given the stress of the previous time step in the material model inputs in, the elastic shear moduli average_elastic_shear_moduli a each point, and the (viscous) viscosities given in the material model outputs object out, compute an update to the elastic stresses and use it to fill the reaction terms material model output property.

§ get_elastic_shear_moduli()

template<int dim>
const std::vector<double>& aspect::MaterialModel::Rheology::Elasticity< dim >::get_elastic_shear_moduli ( ) const

Return the values of the elastic shear moduli for each composition used in the rheology model.

§ calculate_elastic_viscosity()

template<int dim>
double aspect::MaterialModel::Rheology::Elasticity< dim >::calculate_elastic_viscosity ( const double  shear_modulus) const

Calculates the effective elastic viscosity (this is the equivalent viscosity of a material which was unstressed at the end of the previous timestep).

§ calculate_viscoelastic_viscosity()

template<int dim>
double aspect::MaterialModel::Rheology::Elasticity< dim >::calculate_viscoelastic_viscosity ( const double  viscosity,
const double  shear_modulus 
) const

Given the (viscous or visco-plastic) viscosity and the shear modulus, compute the viscoelastic viscosity (eqn 28 in Moresi et al., 2003, J. Comp. Phys.).

§ calculate_viscoelastic_strain_rate()

template<int dim>
SymmetricTensor<2,dim> aspect::MaterialModel::Rheology::Elasticity< dim >::calculate_viscoelastic_strain_rate ( const SymmetricTensor< 2, dim > &  strain_rate,
const SymmetricTensor< 2, dim > &  stored_stress,
const double  shear_modulus 
) const

Calculate the effective deviatoric strain rate tensor, which equals the true deviatoric strain rate plus a fictional strain rate which would arise from stored elastic stresses. In ASPECT, this additional strain rate is supported by a fictional body force. This formulation allows the use of an isotropic effective viscosity by ensuring that the resulting strain rate tensor is equal to the total current stress tensor multiplied by a scalar.

§ elastic_timestep()

template<int dim>
double aspect::MaterialModel::Rheology::Elasticity< dim >::elastic_timestep ( ) const

Compute the elastic time step.

Member Data Documentation

§ elastic_damper_viscosity

template<int dim>
double aspect::MaterialModel::Rheology::Elasticity< dim >::elastic_damper_viscosity
private

Viscosity of a damper used to stabilize elasticity. A value of 0 Pas is equivalent to not using a damper.

Definition at line 156 of file elasticity.h.

§ elastic_shear_moduli

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::Elasticity< dim >::elastic_shear_moduli
private

Vector for field elastic shear moduli, read from parameter file.

Definition at line 161 of file elasticity.h.

§ use_fixed_elastic_time_step

template<int dim>
bool aspect::MaterialModel::Rheology::Elasticity< dim >::use_fixed_elastic_time_step
private

Bool indicating whether to use a fixed material time scale in the viscoelastic rheology for all time steps (if true) or to use the actual (variable) advection time step of the model (if false). Read from parameter file.

Definition at line 169 of file elasticity.h.

§ fixed_elastic_time_step

template<int dim>
double aspect::MaterialModel::Rheology::Elasticity< dim >::fixed_elastic_time_step
private

Double for fixed elastic time step value, read from parameter file

Definition at line 174 of file elasticity.h.

§ stabilization_time_scale_factor

template<int dim>
double aspect::MaterialModel::Rheology::Elasticity< dim >::stabilization_time_scale_factor
private

A stabilization factor for the elastic stresses that influences how fast elastic stresses adjust to deformation. 1.0 is equivalent to no stabilization, and infinity is equivalent to not applying elastic stresses at all. The factor is multiplied with the computational time step to create a time scale.

Definition at line 183 of file elasticity.h.

§ evaluator

template<int dim>
std::unique_ptr<FEPointEvaluation<dim, dim> > aspect::MaterialModel::Rheology::Elasticity< dim >::evaluator
mutableprivate

We cache the evaluator that is necessary to evaluate the old velocity gradients. They are required to compute the elastic stresses, but are not provided by the material model. By caching the evaluator, we can avoid recreating it every time we need it.

Definition at line 192 of file elasticity.h.


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