ASPECT
peierls_creep.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 2020 - 2023 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_peierls_creep_h
22 #define _aspect_material_model_rheology_peierls_creep_h
23 
24 #include <aspect/global.h>
26 
27 namespace aspect
28 {
29  namespace MaterialModel
30  {
31  using namespace dealii;
32 
33  namespace Rheology
34  {
39  {
45  double prefactor;
53  double stress_cutoff;
54 
59  };
60 
69  template <int dim>
71  {
72  public:
76  PeierlsCreep();
77 
82  compute_creep_parameters (const unsigned int composition,
83  const std::vector<double> &phase_function_values = std::vector<double>(),
84  const std::vector<unsigned int> &n_phases_per_composition = std::vector<unsigned int>()) const;
85 
89  static
90  void
91  declare_parameters (ParameterHandler &prm);
92 
100  void
101  parse_parameters (ParameterHandler &prm,
102  const std::unique_ptr<std::vector<unsigned int>> &expected_n_phases_per_composition = nullptr);
103 
111  double
112  compute_approximate_viscosity (const double strain_rate,
113  const double pressure,
114  const double temperature,
115  const unsigned int composition,
116  const std::vector<double> &phase_function_values = std::vector<double>(),
117  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
118 
126  double
127  compute_exact_viscosity (const double strain_rate,
128  const double pressure,
129  const double temperature,
130  const unsigned int composition,
131  const std::vector<double> &phase_function_values = std::vector<double>(),
132  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
133 
143  double
144  compute_viscosity (const double strain_rate,
145  const double pressure,
146  const double temperature,
147  const unsigned int composition,
148  const std::vector<double> &phase_function_values = std::vector<double>(),
149  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
150 
155  std::pair<double, double>
156  compute_approximate_strain_rate_and_derivative (const double stress,
157  const double pressure,
158  const double temperature,
159  const PeierlsCreepParameters creep_parameters) const;
160 
165  std::pair<double, double>
166  compute_exact_strain_rate_and_derivative (const double stress,
167  const double pressure,
168  const double temperature,
169  const PeierlsCreepParameters creep_parameters) const;
170 
176  std::pair<double, double>
177  compute_exact_log_strain_rate_and_derivative (const double log_stress,
178  const double pressure,
179  const double temperature,
180  const PeierlsCreepParameters creep_parameters) const;
181 
187  std::pair<double, double>
188  compute_approximate_log_strain_rate_and_derivative (const double log_stress,
189  const double pressure,
190  const double temperature,
191  const PeierlsCreepParameters creep_parameters) const;
192 
200  std::pair<double, double>
201  compute_strain_rate_and_derivative (const double stress,
202  const double pressure,
203  const double temperature,
204  const PeierlsCreepParameters creep_parameters) const;
205 
206  private:
220  {
222  exact
223  } peierls_creep_flow_law;
224 
228  std::vector<double> prefactors;
229 
233  std::vector<double> stress_exponents;
234 
238  std::vector<double> activation_energies;
239 
243  std::vector<double> activation_volumes;
244 
248  std::vector<double> peierls_stresses;
249 
253  std::vector<double> fitting_parameters;
254 
259  std::vector<double> glide_parameters_p;
260 
265  std::vector<double> glide_parameters_q;
266 
267  std::vector<double> stress_cutoffs;
268 
274 
281 
282  };
283  }
284  }
285 }
286 #endif
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:59
Definition: compat.h:42