ASPECT
peierls_creep.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_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  namespace Rheology
32  {
37  {
43  double prefactor;
51  double stress_cutoff;
52 
57  };
58 
67  template <int dim>
69  {
70  public:
74  PeierlsCreep();
75 
80  compute_creep_parameters (const unsigned int composition,
81  const std::vector<double> &phase_function_values = std::vector<double>(),
82  const std::vector<unsigned int> &n_phases_per_composition = std::vector<unsigned int>()) const;
83 
87  static
88  void
89  declare_parameters (ParameterHandler &prm);
90 
98  void
99  parse_parameters (ParameterHandler &prm,
100  const std::unique_ptr<std::vector<unsigned int>> &expected_n_phases_per_composition = nullptr);
101 
109  double
110  compute_approximate_viscosity (const double strain_rate,
111  const double pressure,
112  const double temperature,
113  const unsigned int composition,
114  const std::vector<double> &phase_function_values = std::vector<double>(),
115  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
116 
124  double
125  compute_exact_viscosity (const double strain_rate,
126  const double pressure,
127  const double temperature,
128  const unsigned int composition,
129  const std::vector<double> &phase_function_values = std::vector<double>(),
130  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
131 
141  double
142  compute_viscosity (const double strain_rate,
143  const double pressure,
144  const double temperature,
145  const unsigned int composition,
146  const std::vector<double> &phase_function_values = std::vector<double>(),
147  const std::vector<unsigned int> &n_phase_transitions_per_composition = std::vector<unsigned int>()) const;
148 
153  std::pair<double, double>
154  compute_approximate_strain_rate_and_derivative (const double stress,
155  const double pressure,
156  const double temperature,
157  const PeierlsCreepParameters creep_parameters) const;
158 
163  std::pair<double, double>
164  compute_exact_strain_rate_and_derivative (const double stress,
165  const double pressure,
166  const double temperature,
167  const PeierlsCreepParameters creep_parameters) const;
168 
174  std::pair<double, double>
175  compute_exact_log_strain_rate_and_derivative (const double log_stress,
176  const double pressure,
177  const double temperature,
178  const PeierlsCreepParameters creep_parameters) const;
179 
185  std::pair<double, double>
186  compute_approximate_log_strain_rate_and_derivative (const double log_stress,
187  const double pressure,
188  const double temperature,
189  const PeierlsCreepParameters creep_parameters) const;
190 
198  std::pair<double, double>
199  compute_strain_rate_and_derivative (const double stress,
200  const double pressure,
201  const double temperature,
202  const PeierlsCreepParameters creep_parameters) const;
203 
204  private:
218  {
220  exact
221  } peierls_creep_flow_law;
222 
226  std::vector<double> prefactors;
227 
231  std::vector<double> stress_exponents;
232 
236  std::vector<double> activation_energies;
237 
241  std::vector<double> activation_volumes;
242 
246  std::vector<double> peierls_stresses;
247 
251  std::vector<double> fitting_parameters;
252 
257  std::vector<double> glide_parameters_p;
258 
263  std::vector<double> glide_parameters_q;
264 
265  std::vector<double> stress_cutoffs;
266 
272 
279 
280  };
281  }
282  }
283 }
284 #endif
void declare_parameters(ParameterHandler &prm)
Definition: compat.h:59