ASPECT
|
Public Member Functions | |
PeierlsCreep () | |
const PeierlsCreepParameters | compute_creep_parameters (const unsigned int composition, const std::vector< double > &phase_function_values=std::vector< double >(), const std::vector< unsigned int > &n_phases_per_composition=std::vector< unsigned int >()) const |
void | parse_parameters (ParameterHandler &prm, const std::unique_ptr< std::vector< unsigned int >> &expected_n_phases_per_composition=nullptr) |
double | compute_approximate_viscosity (const double strain_rate, const double pressure, const double temperature, 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_exact_viscosity (const double strain_rate, const double pressure, const double temperature, 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_viscosity (const double strain_rate, const double pressure, const double temperature, 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 |
std::pair< double, double > | compute_approximate_strain_rate_and_derivative (const double stress, const double pressure, const double temperature, const PeierlsCreepParameters creep_parameters) const |
std::pair< double, double > | compute_exact_strain_rate_and_derivative (const double stress, const double pressure, const double temperature, const PeierlsCreepParameters creep_parameters) const |
std::pair< double, double > | compute_exact_log_strain_rate_and_derivative (const double log_stress, const double pressure, const double temperature, const PeierlsCreepParameters creep_parameters) const |
std::pair< double, double > | compute_approximate_log_strain_rate_and_derivative (const double log_stress, const double pressure, const double temperature, const PeierlsCreepParameters creep_parameters) const |
std::pair< double, double > | compute_strain_rate_and_derivative (const double stress, const double pressure, const double temperature, const PeierlsCreepParameters creep_parameters) 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::BlockVector & | get_current_linearization_point () const |
const LinearAlgebra::BlockVector & | get_solution () const |
const LinearAlgebra::BlockVector & | get_old_solution () const |
const LinearAlgebra::BlockVector & | get_old_old_solution () const |
const LinearAlgebra::BlockVector & | get_reaction_vector () const |
const LinearAlgebra::BlockVector & | get_mesh_velocity () const |
const DoFHandler< dim > & | get_dof_handler () const |
const FiniteElement< dim > & | get_fe () const |
const LinearAlgebra::BlockSparseMatrix & | get_system_matrix () const |
const LinearAlgebra::BlockSparseMatrix & | get_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_managers () const |
const Particle::Manager< dim > & | get_particle_manager (const unsigned int particle_manager_index) const |
Particle::Manager< dim > & | get_particle_manager (const unsigned int particle_manager_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 Types | |
enum | PeierlsCreepScheme { viscosity_approximation, exact } |
Private Attributes | |
enum aspect::MaterialModel::Rheology::PeierlsCreep::PeierlsCreepScheme | peierls_creep_flow_law |
std::vector< double > | prefactors |
std::vector< double > | stress_exponents |
std::vector< double > | activation_energies |
std::vector< double > | activation_volumes |
std::vector< double > | peierls_stresses |
std::vector< double > | fitting_parameters |
std::vector< double > | glide_parameters_p |
std::vector< double > | glide_parameters_q |
std::vector< double > | stress_cutoffs |
bool | apply_strict_cutoff |
double | strain_rate_residual_threshold |
unsigned int | stress_max_iteration_number |
Peierls creep is a low temperature, high stress viscous deformation mechanism. Crystal deformation occurs through dislocation glide, which requires high stresses. It is considered an important deformation mechanism in the lithosphere and subducting slabs (Kumamoto et al., 2017, Science Advances). The approximate form of the Peierls flow law used here is based on a derivation from Kameyama et al., 1999, Earth and Planetary Science Letters.
Definition at line 70 of file peierls_creep.h.
|
private |
Enumeration for selecting which type of Peierls creep flow law to use. The options are currently "exact" (which requires a Newton-Raphson iteration in compute_viscosity) and "viscosity_approximation", which uses a simplified expression where the strain rate has a power law dependence on the stress (so that the equation can be directly inverted for viscosity). This approximation requires specifying one fitting parameter (gamma) to obtain the best fit in the expected stress range for a given problem (gamma = stress / peierls_stress). The derivation for this approximation (derived by Magali Billen) can be found at: https://ucdavis.app.box.com/s/cl5mwhkjeabol4otrdukfcdwfvg9me4w/file/705438695737.
Enumerator | |
---|---|
viscosity_approximation | |
exact |
Definition at line 219 of file peierls_creep.h.
aspect::MaterialModel::Rheology::PeierlsCreep< dim >::PeierlsCreep | ( | ) |
Constructor.
const PeierlsCreepParameters aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_creep_parameters | ( | const unsigned int | composition, |
const std::vector< double > & | phase_function_values = std::vector< double >() , |
||
const std::vector< unsigned int > & | n_phases_per_composition = std::vector< unsigned int >() |
||
) | const |
Compute the creep parameters for the Peierls creep law.
|
static |
Declare the parameters this function takes through input files.
void aspect::MaterialModel::Rheology::PeierlsCreep< dim >::parse_parameters | ( | ParameterHandler & | prm, |
const std::unique_ptr< std::vector< unsigned int >> & | expected_n_phases_per_composition = nullptr |
||
) |
Read the parameters 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 and will be checked against the parsed parameters.
double aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_approximate_viscosity | ( | const double | strain_rate, |
const double | pressure, | ||
const double | temperature, | ||
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 viscosity based on the approximate Peierls creep flow law. 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.
double aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_exact_viscosity | ( | const double | strain_rate, |
const double | pressure, | ||
const double | temperature, | ||
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 viscosity based on the exact Peierls creep flow law. 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.
double aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_viscosity | ( | const double | strain_rate, |
const double | pressure, | ||
const double | temperature, | ||
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 viscosity based on the selected Peierls creep flow law. This function uses either the compute_approximate_viscosity or the compute_exact_viscosity function. 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.
std::pair<double, double> aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_approximate_strain_rate_and_derivative | ( | const double | stress, |
const double | pressure, | ||
const double | temperature, | ||
const PeierlsCreepParameters | creep_parameters | ||
) | const |
Compute the strain rate and first stress derivative as a function of stress based on the approximate Peierls creep law.
std::pair<double, double> aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_exact_strain_rate_and_derivative | ( | const double | stress, |
const double | pressure, | ||
const double | temperature, | ||
const PeierlsCreepParameters | creep_parameters | ||
) | const |
Compute the strain rate and first stress derivative as a function of stress based on the exact Peierls creep law.
std::pair<double, double> aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_exact_log_strain_rate_and_derivative | ( | const double | log_stress, |
const double | pressure, | ||
const double | temperature, | ||
const PeierlsCreepParameters | creep_parameters | ||
) | const |
Compute the natural logarithm of the strain rate norm and its first derivative with respect to the natural logarithm of the stress norm based on the exact Peierls creep law.
std::pair<double, double> aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_approximate_log_strain_rate_and_derivative | ( | const double | log_stress, |
const double | pressure, | ||
const double | temperature, | ||
const PeierlsCreepParameters | creep_parameters | ||
) | const |
Compute the natural logarithm of the strain rate norm and its first derivative with respect to the natural logarithm of the stress norm based on the approximate Peierls creep law.
std::pair<double, double> aspect::MaterialModel::Rheology::PeierlsCreep< dim >::compute_strain_rate_and_derivative | ( | const double | stress, |
const double | pressure, | ||
const double | temperature, | ||
const PeierlsCreepParameters | creep_parameters | ||
) | const |
Compute the strain rate and first stress derivative as a function of stress based on the selected Peierls creep law. This function uses either the compute_approximate_strain_rate_and_derivative or the compute_exact_strain_rate_and_derivative function.
|
private |
|
private |
List of Peierls creep prefactors (A).
Definition at line 228 of file peierls_creep.h.
|
private |
List of Peierls creep stress exponents (n).
Definition at line 233 of file peierls_creep.h.
|
private |
List of Peierls creep activation energies (E).
Definition at line 238 of file peierls_creep.h.
|
private |
List of Peierls creep activation volumes (V).
Definition at line 243 of file peierls_creep.h.
|
private |
List of Peierls stresses (sigma_p).
Definition at line 248 of file peierls_creep.h.
|
private |
List of Peierls fitting parameters (gamma).
Definition at line 253 of file peierls_creep.h.
|
private |
List of the first Peierls parameter related to dislocation glide (p).
Definition at line 259 of file peierls_creep.h.
|
private |
List of the second Peierls parameter related to dislocation glide (q).
Definition at line 265 of file peierls_creep.h.
|
private |
Definition at line 267 of file peierls_creep.h.
|
private |
A parameter determines whether a strict cutoff on the stress is applied to the Peierls creep
Definition at line 273 of file peierls_creep.h.
|
private |
Parameters governing the iteration for the exact Peierls viscosity.
Definition at line 279 of file peierls_creep.h.
|
private |
Definition at line 280 of file peierls_creep.h.