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

Public Member Functions

 DruckerPragerPower ()
 
void parse_parameters (ParameterHandler &prm, const std::unique_ptr< std::vector< unsigned int >> &expected_n_phases_per_composition=nullptr)
 
const DruckerPragerParameters compute_drucker_prager_parameters (const unsigned int composition, const std::vector< double > &phase_function_values=std::vector< double >(), const std::vector< unsigned int > &n_phase_transitions_per_composition=std::vector< unsigned int >()) const
 
double compute_yield_stress (const double cohesion, const double angle_internal_friction, const double pressure, const double max_yield_stress) const
 
double compute_viscosity (const double cohesion, const double angle_internal_friction, const double pressure, const double effective_strain_rate, const double max_yield_stress) const
 
std::pair< double, double > compute_strain_rate_and_derivative (const double stress, const double pressure, const DruckerPragerParameters &p) const
 
std::pair< double, double > compute_log_strain_rate_and_derivative (const double log_stress, const double pressure, const DruckerPragerParameters &p) 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)
 

Private Attributes

std::vector< double > angles_internal_friction
 
std::vector< double > cohesions
 
double max_yield_stress
 
double drucker_prager_edot_ref
 
double drucker_prager_log_edot_ref
 
double drucker_prager_stress_exponent
 

Detailed Description

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

Definition at line 41 of file drucker_prager_power.h.

Constructor & Destructor Documentation

§ DruckerPragerPower()

The Drucker-Prager power-law rheology approximates the standard Drucker-Prager plastic model. The standard model involves a simple plastic element that yields at a stress of (6.0 * cohesion * cos_phi + 6.0 * pressure * sin_phi) / (sqrt(3.0) * (3.0 + sin_phi)) in 3D or cohesion * cos_phi + pressure * sin_phi in 2D. The "Power-Law" modification introduced in this class replaces the yield criterion with a power law dependence of stress on the strain rate: stress = yield_stress * (edot/edot_ref)^(1/n). In the limit that n -> infinity, the model approaches the classic Drucker-Prager model.

Member Function Documentation

§ declare_parameters()

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

Declare the parameters this function takes through input files.

§ parse_parameters()

template<int dim>
void aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::parse_parameters ( ParameterHandler &  prm,
const std::unique_ptr< std::vector< unsigned int >> &  expected_n_phases_per_composition = nullptr 
)

Read the parameters this class declares from the parameter file. If expected_n_phases_per_composition points to a vector of unsigned integers, this is considered the number of phases for each compositional field (plus possibly a background field) and this number will be checked against the parsed parameters.

Parameters
[in]prmThe ParameterHandler to read from.
expected_n_phases_per_compositionOptional list of number of phases.

§ compute_drucker_prager_parameters()

template<int dim>
const DruckerPragerParameters aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::compute_drucker_prager_parameters ( const unsigned int  composition,
const std::vector< double > &  phase_function_values = std::vector< double >(),
const std::vector< unsigned int > &  n_phase_transitions_per_composition = std::vector< unsigned int >() 
) const

Compute the parameters for the Drucker Prager plasticity. If n_phase_transitions_per_composition points to a vector of unsigned integers this is considered the number of phase transitions for each compositional field and viscosity will be first computed on each phase and then averaged for each compositional field.

§ compute_yield_stress()

template<int dim>
double aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::compute_yield_stress ( const double  cohesion,
const double  angle_internal_friction,
const double  pressure,
const double  max_yield_stress 
) const

Compute the plastic yield stress based on the Drucker Prager yield criterion.

§ compute_viscosity()

template<int dim>
double aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::compute_viscosity ( const double  cohesion,
const double  angle_internal_friction,
const double  pressure,
const double  effective_strain_rate,
const double  max_yield_stress 
) const

Compute the apparent viscosity using the yield stress and effective strain rate.

§ compute_strain_rate_and_derivative()

template<int dim>
std::pair<double, double> aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::compute_strain_rate_and_derivative ( const double  stress,
const double  pressure,
const DruckerPragerParameters p 
) const

Compute the strain rate and first derivative as a function of stress according to the damped Drucker-Prager power-law flow law.

§ compute_log_strain_rate_and_derivative()

template<int dim>
std::pair<double, double> aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::compute_log_strain_rate_and_derivative ( const double  log_stress,
const double  pressure,
const DruckerPragerParameters p 
) const

Compute the natural logarithm of the strain rate and first log stress derivative as a function of log stress according to the damped Drucker-Prager power-law flow law.

Member Data Documentation

§ angles_internal_friction

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::angles_internal_friction
private

The angles of internal friction (phi) are input by the user in degrees, but stored as radians.

Definition at line 134 of file drucker_prager_power.h.

§ cohesions

template<int dim>
std::vector<double> aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::cohesions
private

The cohesion is provided and stored in Pa.

Definition at line 139 of file drucker_prager_power.h.

§ max_yield_stress

template<int dim>
double aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::max_yield_stress
private

The yield stress is limited to a constant value, stored in Pa.

Definition at line 144 of file drucker_prager_power.h.

§ drucker_prager_edot_ref

template<int dim>
double aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::drucker_prager_edot_ref
private

The reference strain rate at which the stress is equal to the "yield stress".

Definition at line 150 of file drucker_prager_power.h.

§ drucker_prager_log_edot_ref

template<int dim>
double aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::drucker_prager_log_edot_ref
private

The natural logarithm of the reference strain rate, stored for convenience.

Definition at line 156 of file drucker_prager_power.h.

§ drucker_prager_stress_exponent

template<int dim>
double aspect::MaterialModel::Rheology::DruckerPragerPower< dim >::drucker_prager_stress_exponent
private

The stress exponent n of the pseudo-plastic element.

Definition at line 161 of file drucker_prager_power.h.


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