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  using namespace dealii;
38 
39  namespace Rheology
40  {
49  namespace ViscosityAveraging
50  {
51  enum Kind
52  {
55  };
56 
63  parse (const std::string &parameter_name,
64  const ParameterHandler &prm);
65  }
66 
67  template <int dim>
69  {
70  public:
75 
79  static
80  void
81  declare_parameters (ParameterHandler &prm);
82 
90  void
91  parse_parameters (ParameterHandler &prm,
92  const std::unique_ptr<std::vector<unsigned int>> &expected_n_phases_per_composition = nullptr);
93 
101  double
102  compute_viscosity (const double pressure,
103  const double temperature,
104  const double grain_size,
105  const std::vector<double> &volume_fractions,
106  const SymmetricTensor<2,dim> &strain_rate,
107  std::vector<double> &partial_strain_rates,
108  const std::vector<double> &phase_function_values = std::vector<double>(),
109  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
110 
111  private:
120  double
121  compute_isostress_viscosity (const double pressure,
122  const double temperature,
123  const double grain_size,
124  const std::vector<double> &volume_fractions,
125  const SymmetricTensor<2,dim> &strain_rate,
126  std::vector<double> &partial_strain_rates,
127  const std::vector<double> &phase_function_values = std::vector<double>(),
128  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
129 
135  std::pair<double, double>
136  calculate_isostress_log_strain_rate_and_derivative(const std::vector<std::array<std::pair<double, double>, 4>> &logarithmic_strain_rates_and_stress_derivatives,
137  const double viscoplastic_stress,
138  std::vector<double> &partial_strain_rates) const;
139 
148  double
149  compute_isostrain_viscosity (const double pressure,
150  const double temperature,
151  const double grain_size,
152  const std::vector<double> &volume_fractions,
153  const SymmetricTensor<2,dim> &strain_rate,
154  std::vector<double> &partial_strain_rates,
155  const std::vector<double> &phase_function_values = std::vector<double>(),
156  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
157 
158 
167  double
168  compute_composition_viscosity (const double pressure,
169  const double temperature,
170  const double grain_size,
171  const unsigned int composition,
172  const SymmetricTensor<2,dim> &strain_rate,
173  std::vector<double> &partial_strain_rates,
174  const std::vector<double> &phase_function_values = std::vector<double>(),
175  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
176 
183  std::pair<double, double>
184  calculate_composition_log_strain_rate_and_derivative(const std::array<std::pair<double, double>, 4> &logarithmic_strain_rates_and_stress_derivatives,
185  const double viscoplastic_stress,
186  std::vector<double> &partial_strain_rates) const;
187 
192 
200 
204  std::vector<unsigned int> active_flow_mechanisms;
205 
212  static constexpr unsigned int n_decomposed_strain_rates = 5;
213 
218  std::unique_ptr<Rheology::DiffusionCreep<dim>> diffusion_creep;
219  std::unique_ptr<Rheology::DislocationCreep<dim>> dislocation_creep;
220  std::unique_ptr<Rheology::PeierlsCreep<dim>> peierls_creep;
221  std::unique_ptr<Rheology::DruckerPragerPower<dim>> drucker_prager;
222 
227 
233 
234 
244 
253 
258 
265 
271 
275  const double logmin = std::log(std::numeric_limits<double>::min());
276  };
277  }
278  }
279 }
280 #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
Definition: compat.h:42