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

Public Member Functions

void evaluate (const MaterialModel::MaterialModelInputs< dim > &in, MaterialModel::MaterialModelOutputs< dim > &out) const override
 
Qualitative properties one can ask a material model
bool is_compressible () const override
 
Reference quantities
double reference_viscosity () const override
 
- Public Member Functions inherited from aspect::MaterialModel::Interface< dim >
virtual ~Interface ()
 
virtual void initialize ()
 
virtual void update ()
 
virtual void create_additional_named_outputs (MaterialModelOutputs &outputs) const
 
virtual void fill_additional_material_model_inputs (MaterialModel::MaterialModelInputs< dim > &input, const LinearAlgebra::BlockVector &solution, const FEValuesBase< dim > &fe_values, const Introspection< dim > &introspection) const
 
const NonlinearDependence::ModelDependenceget_model_dependence () const
 
- Public Member Functions inherited from aspect::SimulatorAccess< dim >
 SimulatorAccess ()
 
 SimulatorAccess (const Simulator< dim > &simulator_object)
 
virtual ~SimulatorAccess ()
 
virtual void initialize_simulator (const Simulator< dim > &simulator_object)
 
template<typename PostprocessorType >
PostprocessorType * find_postprocessor () const
 
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
 
TimerOutputget_computing_timer () const
 
const ConditionalOStreamget_pcout () const
 
double get_time () const
 
double get_timestep () const
 
double get_old_timestep () const
 
unsigned int get_timestep_number () 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
 
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
 
void compute_material_model_input_values (const LinearAlgebra::BlockVector &input_solution, const FEValuesBase< dim, dim > &input_finite_element_values, const typename DoFHandler< dim >::active_cell_iterator &cell, const bool compute_strainrate, MaterialModel::MaterialModelInputs< dim > &material_model_inputs) 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
 
DEAL_II_DEPRECATED const BoundaryTemperature::Interface< dim > & get_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
 
DEAL_II_DEPRECATED const BoundaryComposition::Interface< dim > & get_boundary_composition () const
 
const BoundaryComposition::Manager< dim > & get_boundary_composition_manager () const
 
const std::map< types::boundary_id, std::unique_ptr< BoundaryTraction::Interface< dim > > > & get_boundary_traction () const
 
DEAL_II_DEPRECATED const InitialTemperature::Interface< dim > & get_initial_temperature () const
 
const InitialTemperature::Manager< dim > & get_initial_temperature_manager () const
 
DEAL_II_DEPRECATED const InitialComposition::Interface< dim > & get_initial_composition () 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
 
TableHandlerget_statistics_object () const
 
template<typename PostprocessorType >
DEAL_II_DEPRECATED PostprocessorType * find_postprocessor () 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
 

Private Attributes

double reference_T
 
double eta
 
double composition_viscosity_prefactor
 
double thermal_viscosity_exponent
 
double maximum_thermal_prefactor
 
double minimum_thermal_prefactor
 
double k_value
 
EquationOfState::LinearizedIncompressible< dim > equation_of_state
 

Functions used in dealing with run-time parameters

void parse_parameters (ParameterHandler &prm) override
 
static void declare_parameters (ParameterHandler &prm)
 

Additional Inherited Members

- Public Types inherited from aspect::MaterialModel::Interface< dim >
using MaterialModelInputs = MaterialModel::MaterialModelInputs< dim >
 
using MaterialModelOutputs = MaterialModel::MaterialModelOutputs< dim >
 
- Static Public Member Functions inherited from aspect::MaterialModel::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)
 
- Protected Attributes inherited from aspect::MaterialModel::Interface< dim >
NonlinearDependence::ModelDependence model_dependence
 

Detailed Description

template<int dim>
class aspect::MaterialModel::Simple< dim >

A material model that consists of globally constant values for all material parameters except density and viscosity.

The model is considered incompressible, following the definition described in Interface::is_compressible. This is essentially the material model used in the step-32 tutorial program.

Definition at line 45 of file simple.h.

Member Function Documentation

§ evaluate()

template<int dim>
void aspect::MaterialModel::Simple< dim >::evaluate ( const MaterialModel::MaterialModelInputs< dim > &  in,
MaterialModel::MaterialModelOutputs< dim > &  out 
) const
overridevirtual

Function to compute the material properties in out given the inputs in in. If MaterialModelInputs.strain_rate has the length 0, then the viscosity does not need to be computed.

Implements aspect::MaterialModel::Interface< dim >.

§ is_compressible()

template<int dim>
bool aspect::MaterialModel::Simple< dim >::is_compressible ( ) const
overridevirtual

Return whether the model is compressible or not. Incompressibility does not necessarily imply that the density is constant; rather, it may still depend on temperature or pressure. In the current context, compressibility means whether we should solve the continuity equation as \(\nabla \cdot (\rho \mathbf u)=0\) (compressible Stokes) or as \(\nabla \cdot \mathbf{u}=0\) (incompressible Stokes).

Implements aspect::MaterialModel::Interface< dim >.

§ reference_viscosity()

template<int dim>
double aspect::MaterialModel::Simple< dim >::reference_viscosity ( ) const
overridevirtual

Return a reference value typical of the viscosities that appear in this model. This value is not actually used in the material description itself, but is used in scaling variables to the same numerical order of magnitude when solving linear systems. Specifically, the reference viscosity appears in the factor scaling the pressure against the velocity. It is also used in computing dimension-less quantities. You may want to take a look at the Kronbichler, Heister, Bangerth 2012 paper that describes the design of ASPECT for a description of this pressure scaling.

Note
The reference viscosity should take into account the complete constitutive relationship, defined as the scalar viscosity times the constitutive tensor. In most cases, the constitutive tensor will simply be the identity tensor (this is the default case), but this may become important for material models with anisotropic viscosities, if the constitutive tensor is not normalized.

Implements aspect::MaterialModel::Interface< dim >.

§ declare_parameters()

template<int dim>
static void aspect::MaterialModel::Simple< dim >::declare_parameters ( ParameterHandler prm)
static

Declare the parameters this class takes through input files.

§ parse_parameters()

template<int dim>
void aspect::MaterialModel::Simple< dim >::parse_parameters ( ParameterHandler prm)
overridevirtual

Read the parameters this class declares from the parameter file.

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

Member Data Documentation

§ reference_T

template<int dim>
double aspect::MaterialModel::Simple< dim >::reference_T
private

Definition at line 100 of file simple.h.

§ eta

template<int dim>
double aspect::MaterialModel::Simple< dim >::eta
private

Definition at line 101 of file simple.h.

§ composition_viscosity_prefactor

template<int dim>
double aspect::MaterialModel::Simple< dim >::composition_viscosity_prefactor
private

Definition at line 102 of file simple.h.

§ thermal_viscosity_exponent

template<int dim>
double aspect::MaterialModel::Simple< dim >::thermal_viscosity_exponent
private

Definition at line 103 of file simple.h.

§ maximum_thermal_prefactor

template<int dim>
double aspect::MaterialModel::Simple< dim >::maximum_thermal_prefactor
private

Definition at line 104 of file simple.h.

§ minimum_thermal_prefactor

template<int dim>
double aspect::MaterialModel::Simple< dim >::minimum_thermal_prefactor
private

Definition at line 105 of file simple.h.

§ k_value

template<int dim>
double aspect::MaterialModel::Simple< dim >::k_value
private

The thermal conductivity.

Definition at line 110 of file simple.h.

§ equation_of_state

template<int dim>
EquationOfState::LinearizedIncompressible<dim> aspect::MaterialModel::Simple< dim >::equation_of_state
private

Definition at line 112 of file simple.h.


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