21 #ifndef _aspect_material_model_rheology_composite_visco_plastic_h 22 #define _aspect_material_model_rheology_composite_visco_plastic_h 35 namespace MaterialModel
47 namespace ViscosityAveraging
61 parse (
const std::string ¶meter_name,
62 const ParameterHandler &prm);
89 parse_parameters (ParameterHandler &prm,
90 const std::unique_ptr<std::vector<unsigned int>> &expected_n_phases_per_composition =
nullptr);
100 compute_viscosity (
const double pressure,
102 const double grain_size,
103 const std::vector<double> &volume_fractions,
105 std::vector<double> &partial_strain_rates,
106 const std::vector<double> &phase_function_values = std::vector<double>(),
107 const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>())
const;
119 compute_isostress_viscosity (
const double pressure,
120 const double temperature,
121 const double grain_size,
122 const std::vector<double> &volume_fractions,
123 const SymmetricTensor<2,dim> &strain_rate,
124 std::vector<double> &partial_strain_rates,
125 const std::vector<double> &phase_function_values = std::vector<double>(),
126 const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>())
const;
133 std::pair<double, double>
134 calculate_isostress_log_strain_rate_and_derivative(
const std::vector<std::array<std::pair<double, double>, 4>> &logarithmic_strain_rates_and_stress_derivatives,
135 const double viscoplastic_stress,
136 std::vector<double> &partial_strain_rates)
const;
147 compute_isostrain_viscosity (
const double pressure,
148 const double temperature,
149 const double grain_size,
150 const std::vector<double> &volume_fractions,
151 const SymmetricTensor<2,dim> &strain_rate,
152 std::vector<double> &partial_strain_rates,
153 const std::vector<double> &phase_function_values = std::vector<double>(),
154 const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>())
const;
166 compute_composition_viscosity (
const double pressure,
167 const double temperature,
168 const double grain_size,
169 const unsigned int composition,
170 const SymmetricTensor<2,dim> &strain_rate,
171 std::vector<double> &partial_strain_rates,
172 const std::vector<double> &phase_function_values = std::vector<double>(),
173 const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>())
const;
181 std::pair<double, double>
182 calculate_composition_log_strain_rate_and_derivative(
const std::array<std::pair<double, double>, 4> &logarithmic_strain_rates_and_stress_derivatives,
183 const double viscoplastic_stress,
184 std::vector<double> &partial_strain_rates)
const;
210 static constexpr
unsigned int n_decomposed_strain_rates = 5;
273 const double logmin = std::log(std::numeric_limits<double>::min());
std::vector< unsigned int > active_flow_mechanisms
std::unique_ptr< Rheology::PeierlsCreep< dim > > peierls_creep
double log_strain_rate_residual_threshold
double strain_rate_scaling_factor
std::unique_ptr< Rheology::DruckerPragerPower< dim > > drucker_prager
bool use_dislocation_creep
void declare_parameters(ParameterHandler &prm)
std::unique_ptr< Rheology::DislocationCreep< dim > > dislocation_creep
unsigned int number_of_chemical_compositions
ViscosityAveraging::Kind parse(const std::string ¶meter_name, const ParameterHandler &prm)
double minimum_strain_rate
std::unique_ptr< Rheology::DiffusionCreep< dim > > diffusion_creep
unsigned int stress_max_iteration_number
ViscosityAveraging::Kind viscosity_averaging_scheme