ASPECT
composite_visco_plastic.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2020 - 2024 by the authors of the ASPECT code.
3 
4  This file is part of ASPECT.
5 
6  ASPECT is free software; you can redistribute it and/or modify
7  it under the terms of the GNU General Public License as published by
8  the Free Software Foundation; either version 2, or (at your option)
9  any later version.
10 
11  ASPECT is distributed in the hope that it will be useful,
12  but WITHOUT ANY WARRANTY; without even the implied warranty of
13  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
14  GNU General Public License for more details.
15 
16  You should have received a copy of the GNU General Public License
17  along with ASPECT; see the file LICENSE. If not see
18  <http://www.gnu.org/licenses/>.
19 */
20 
21 #ifndef _aspect_material_model_rheology_composite_visco_plastic_h
22 #define _aspect_material_model_rheology_composite_visco_plastic_h
23 
24 #include <aspect/global.h>
32 
33 namespace aspect
34 {
35  namespace MaterialModel
36  {
37  namespace Rheology
38  {
47  namespace ViscosityAveraging
48  {
49  enum Kind
50  {
53  };
54 
61  parse (const std::string &parameter_name,
62  const ParameterHandler &prm);
63  }
64 
65  template <int dim>
67  {
68  public:
73 
77  static
78  void
79  declare_parameters (ParameterHandler &prm);
80 
88  void
89  parse_parameters (ParameterHandler &prm,
90  const std::unique_ptr<std::vector<unsigned int>> &expected_n_phases_per_composition = nullptr);
91 
99  double
100  compute_viscosity (const double pressure,
101  const double temperature,
102  const double grain_size,
103  const std::vector<double> &volume_fractions,
104  const SymmetricTensor<2,dim> &strain_rate,
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;
108 
109  private:
118  double
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;
127 
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;
137 
146  double
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;
155 
156 
165  double
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;
174 
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;
185 
190 
198 
202  std::vector<unsigned int> active_flow_mechanisms;
203 
210  static constexpr unsigned int n_decomposed_strain_rates = 5;
211 
216  std::unique_ptr<Rheology::DiffusionCreep<dim>> diffusion_creep;
217  std::unique_ptr<Rheology::DislocationCreep<dim>> dislocation_creep;
218  std::unique_ptr<Rheology::PeierlsCreep<dim>> peierls_creep;
219  std::unique_ptr<Rheology::DruckerPragerPower<dim>> drucker_prager;
220 
225 
231 
232 
242 
251 
256 
263 
269 
273  const double logmin = std::log(std::numeric_limits<double>::min());
274  };
275  }
276  }
277 }
278 #endif
std::unique_ptr< Rheology::PeierlsCreep< dim > > peierls_creep
std::unique_ptr< Rheology::DruckerPragerPower< dim > > drucker_prager
void declare_parameters(ParameterHandler &prm)
std::unique_ptr< Rheology::DislocationCreep< dim > > dislocation_creep
ViscosityAveraging::Kind parse(const std::string &parameter_name, const ParameterHandler &prm)
std::unique_ptr< Rheology::DiffusionCreep< dim > > diffusion_creep
Definition: compat.h:59