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

Public Member Functions

void parse_parameters (ParameterHandler &prm)
 
std::array< double, 3 > compute_strain_weakening_factors (const std::vector< double > &composition, const unsigned int j) const
 
DEAL_II_DEPRECATED std::array< double, 3 > compute_strain_weakening_factors (const unsigned int j, const std::vector< double > &composition) const
 
std::array< double, 3 > apply_temperature_dependence_to_strain_weakening_factors (const std::array< double, 3 > &weakening_factors, const double temperature, const unsigned int j) const
 
double calculate_strain_healing (const MaterialModel::MaterialModelInputs< dim > &in, const unsigned int j) const
 
std::pair< double, double > calculate_plastic_weakening (const double strain_ii, const unsigned int j) const
 
double calculate_viscous_weakening (const double strain_ii, const unsigned int j) const
 
void compute_finite_strain_reaction_terms (const MaterialModel::MaterialModelInputs< dim > &in, MaterialModel::MaterialModelOutputs< dim > &out) const
 
void fill_reaction_outputs (const MaterialModel::MaterialModelInputs< dim > &in, const int i, const double min_strain_rate, const bool plastic_yielding, MaterialModel::MaterialModelOutputs< dim > &out) const
 
ComponentMask get_strain_composition_mask () const
 
WeakeningMechanism get_weakening_mechanism () const
 
HealingMechanism get_healing_mechanism () 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
 
unsigned int n_particle_worlds () const
 
const Particle::World< dim > & get_particle_world (const unsigned int particle_world_index) const
 
Particle::World< dim > & get_particle_world (const unsigned int particle_world_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 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)
 

Public Attributes

bool use_temperature_activated_strain_softening
 

Private Attributes

WeakeningMechanism weakening_mechanism
 
HealingMechanism healing_mechanism
 
std::vector< double > start_plastic_strain_weakening_intervals
 
std::vector< double > end_plastic_strain_weakening_intervals
 
std::vector< double > cohesion_strain_weakening_factors
 
std::vector< double > friction_strain_weakening_factors
 
std::vector< double > start_viscous_strain_weakening_intervals
 
std::vector< double > end_viscous_strain_weakening_intervals
 
std::vector< double > viscous_strain_weakening_factors
 
std::vector< double > viscous_strain_weakening_T0
 
std::vector< double > viscous_strain_weakening_T1
 
std::vector< double > viscous_strain_weakening_T2
 
std::vector< double > viscous_strain_weakening_T3
 
double strain_healing_temperature_dependent_recovery_rate
 
double strain_healing_temperature_dependent_prefactor
 
std::unique_ptr< FEPointEvaluation< dim, dim > > evaluator
 
std::vector< std::unique_ptr< FEPointEvaluation< 1, dim > > > composition_evaluators
 

Detailed Description

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

Definition at line 74 of file strain_dependent.h.

Member Function Documentation

§ declare_parameters()

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

Declare the parameters this function takes through input files.

§ parse_parameters()

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

Read the parameters from the parameter file.

§ compute_strain_weakening_factors() [1/2]

template<int dim>
std::array<double, 3> aspect::MaterialModel::Rheology::StrainDependent< dim >::compute_strain_weakening_factors ( const std::vector< double > &  composition,
const unsigned int  j 
) const

A function that computes by how much the rheologic parameters change if strain weakening is applied. Given a vector composition of all fields, it returns reduction factors for the cohesion, friction angle and the prefactor of the viscous flow law(s) for the compositional field with index j. The reason all fields are passed is that the weakening factors can depend on the values of fields that track different measures of previously applied strain.

§ compute_strain_weakening_factors() [2/2]

template<int dim>
DEAL_II_DEPRECATED std::array<double, 3> aspect::MaterialModel::Rheology::StrainDependent< dim >::compute_strain_weakening_factors ( const unsigned int  j,
const std::vector< double > &  composition 
) const
Deprecated:
: Deprecated version of the function of the same name described above.

§ apply_temperature_dependence_to_strain_weakening_factors()

template<int dim>
std::array<double, 3> aspect::MaterialModel::Rheology::StrainDependent< dim >::apply_temperature_dependence_to_strain_weakening_factors ( const std::array< double, 3 > &  weakening_factors,
const double  temperature,
const unsigned int  j 
) const

A function that alters the viscous weakening factor based on the temperature field.

§ calculate_strain_healing()

template<int dim>
double aspect::MaterialModel::Rheology::StrainDependent< dim >::calculate_strain_healing ( const MaterialModel::MaterialModelInputs< dim > &  in,
const unsigned int  j 
) const

A function that computes the strain healing (reduction in accumulated strain)

§ calculate_plastic_weakening()

template<int dim>
std::pair<double, double> aspect::MaterialModel::Rheology::StrainDependent< dim >::calculate_plastic_weakening ( const double  strain_ii,
const unsigned int  j 
) const

A function that computes by how much the cohesion and internal friction angle for a given compositional field are weakened under the influence of a given strain.

§ calculate_viscous_weakening()

template<int dim>
double aspect::MaterialModel::Rheology::StrainDependent< dim >::calculate_viscous_weakening ( const double  strain_ii,
const unsigned int  j 
) const

A function that computes by how much the diffusion and dislocation prefactors for a given compositional field are weakened under the influence of a given strain.

§ compute_finite_strain_reaction_terms()

template<int dim>
void aspect::MaterialModel::Rheology::StrainDependent< dim >::compute_finite_strain_reaction_terms ( const MaterialModel::MaterialModelInputs< dim > &  in,
MaterialModel::MaterialModelOutputs< dim > &  out 
) const

A function that fills the reaction terms for the finite strain tensor in MaterialModelOutputs object that is handed over. It assumes the first component of the finite strain tensor is named 's11' and all other components follow this compositional field.

§ fill_reaction_outputs()

template<int dim>
void aspect::MaterialModel::Rheology::StrainDependent< dim >::fill_reaction_outputs ( const MaterialModel::MaterialModelInputs< dim > &  in,
const int  i,
const double  min_strain_rate,
const bool  plastic_yielding,
MaterialModel::MaterialModelOutputs< dim > &  out 
) const

A function that fills the reaction terms for the finite strain invariant(s) in MaterialModelOutputs object that is handed over.

§ get_strain_composition_mask()

template<int dim>
ComponentMask aspect::MaterialModel::Rheology::StrainDependent< dim >::get_strain_composition_mask ( ) const

A function that returns a ComponentMask, which indicates that components associated with strain should be excluded during the volume fraction computation.

§ get_weakening_mechanism()

template<int dim>
WeakeningMechanism aspect::MaterialModel::Rheology::StrainDependent< dim >::get_weakening_mechanism ( ) const

A function that returns the selected type of strain weakening mechanism.

§ get_healing_mechanism()

template<int dim>
HealingMechanism aspect::MaterialModel::Rheology::StrainDependent< dim >::get_healing_mechanism ( ) const

A function that returns the selected type of strain healing mechanism.

Member Data Documentation

§ use_temperature_activated_strain_softening

template<int dim>
bool aspect::MaterialModel::Rheology::StrainDependent< dim >::use_temperature_activated_strain_softening

Whether to use the temperature-activated viscous strain weakening.

Definition at line 149 of file strain_dependent.h.

§ weakening_mechanism

template<int dim>
WeakeningMechanism aspect::MaterialModel::Rheology::StrainDependent< dim >::weakening_mechanism
private

Definition at line 189 of file strain_dependent.h.

§ healing_mechanism

template<int dim>
HealingMechanism aspect::MaterialModel::Rheology::StrainDependent< dim >::healing_mechanism
private

Definition at line 191 of file strain_dependent.h.

§ start_plastic_strain_weakening_intervals

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::start_plastic_strain_weakening_intervals
private

The start of the strain interval (plastic or total strain) within which cohesion and angle of friction should be weakened.

Definition at line 197 of file strain_dependent.h.

§ end_plastic_strain_weakening_intervals

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::end_plastic_strain_weakening_intervals
private

The end of the strain interval (plastic or total strain) within which cohesion and angle of friction should be weakened.

Definition at line 203 of file strain_dependent.h.

§ cohesion_strain_weakening_factors

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::cohesion_strain_weakening_factors
private

The factor specifying the amount of weakening of the cohesion over the prescribed strain interval (plastic or total strain).

Definition at line 209 of file strain_dependent.h.

§ friction_strain_weakening_factors

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::friction_strain_weakening_factors
private

The factor specifying the amount of weakening of the internal friction angles over the prescribed strain interval (plastic or total strain).

Definition at line 216 of file strain_dependent.h.

§ start_viscous_strain_weakening_intervals

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::start_viscous_strain_weakening_intervals
private

The start of the strain interval (viscous or total strain) within which cohesion and angle of friction should be weakened.

Definition at line 222 of file strain_dependent.h.

§ end_viscous_strain_weakening_intervals

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::end_viscous_strain_weakening_intervals
private

The end of the strain interval (viscous or total strain) within which cohesion and angle of friction should be weakened.

Definition at line 228 of file strain_dependent.h.

§ viscous_strain_weakening_factors

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::viscous_strain_weakening_factors
private

The factor specifying the amount of weakening over the prescribed strain interval (viscous or total strain).

Definition at line 234 of file strain_dependent.h.

§ viscous_strain_weakening_T0

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::viscous_strain_weakening_T0
private

The four temperatures that parameterize the temperature-activated strain softening. These can be different for each compositional field. ---— -----— 1 \ / \ / ______/ _ _ _ _ _ viscous_weakening_factor[2]

-------------------------—> T T0 T1 T2 T3

Definition at line 247 of file strain_dependent.h.

§ viscous_strain_weakening_T1

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::viscous_strain_weakening_T1
private

Definition at line 248 of file strain_dependent.h.

§ viscous_strain_weakening_T2

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::viscous_strain_weakening_T2
private

Definition at line 249 of file strain_dependent.h.

§ viscous_strain_weakening_T3

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::StrainDependent< dim >::viscous_strain_weakening_T3
private

Definition at line 250 of file strain_dependent.h.

§ strain_healing_temperature_dependent_recovery_rate

template<int dim>
double aspect::MaterialModel::Rheology::StrainDependent< dim >::strain_healing_temperature_dependent_recovery_rate
private

The healing rate used in the temperature dependent strain healing model.

Definition at line 255 of file strain_dependent.h.

§ strain_healing_temperature_dependent_prefactor

template<int dim>
double aspect::MaterialModel::Rheology::StrainDependent< dim >::strain_healing_temperature_dependent_prefactor
private

A prefactor of viscosity used in the strain healing calculation.

Definition at line 260 of file strain_dependent.h.

§ evaluator

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

We cache the evaluators that are necessary to evaluate the velocity gradients and compositions. By caching the evaluator, we can avoid recreating them every time we need it.

Definition at line 268 of file strain_dependent.h.

§ composition_evaluators

template<int dim>
std::vector<std::unique_ptr<FEPointEvaluation<1, dim> > > aspect::MaterialModel::Rheology::StrainDependent< dim >::composition_evaluators
mutableprivate

Definition at line 269 of file strain_dependent.h.


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