ASPECT
grain_size_evolution.h
Go to the documentation of this file.
1 /*
2  Copyright (C) 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_reaction_model_grain_size_evolution_h
22 #define _aspect_material_model_reaction_model_grain_size_evolution_h
23 
26 
27 namespace aspect
28 {
29  namespace MaterialModel
30  {
31  namespace ReactionModel
32  {
49  template <int dim>
51  {
52  public:
59 
79  void
81  const std::vector<double> &adiabatic_pressure,
82  const std::vector<unsigned int> &phase_indices,
83  const std::function<double(const double,const double,
84  const double,const SymmetricTensor<2,dim> &,const unsigned int,const double, const double)> &dislocation_viscosity,
85  const std::function<double(
86  const double, const double, const double,const double,const double,const unsigned int)> &diffusion_viscosity,
87  const double min_eta,
88  const double max_eta,
89  typename Interface<dim>::MaterialModelOutputs &out) const;
90 
97  void
99 
106  void
108  const typename MaterialModel::MaterialModelOutputs<dim> &out,
109  const std::vector<unsigned int> &phase_indices,
110  const std::vector<double> &dislocation_viscosities,
111  std::vector<std::unique_ptr<MaterialModel::AdditionalMaterialOutputs<dim>>> &additional_outputs) const;
112 
116  static
117  void
118  declare_parameters (ParameterHandler &prm);
119 
123  void
124  parse_parameters (ParameterHandler &prm);
125 
126  private:
130  std::vector<double> grain_growth_activation_energy;
131  std::vector<double> grain_growth_activation_volume;
132  std::vector<double> grain_growth_rate_constant;
133  std::vector<double> grain_growth_exponent;
135  std::vector<double> reciprocal_required_strain;
136  std::vector<double> recrystallized_grain_size;
137 
141  std::vector<double> grain_boundary_energy;
143  std::vector<double> geometric_constant;
144 
149  struct Formulation
150  {
171  enum Kind
172  {
176  };
177 
182  static
183  Kind
184  parse(const std::string &input)
185  {
186  if (input == "paleowattmeter")
188  else if (input == "paleopiezometer")
190  else if (input == "pinned grain damage")
192  else
193  AssertThrow(false, ExcNotImplemented());
194 
195  return Formulation::Kind();
196  }
197  };
198 
204 
211  double compute_partitioning_fraction (const double temperature) const;
212 
222 
234 
244 
254  std::shared_ptr<MaterialUtilities::PhaseFunction<dim>> phase_function;
255 
259  unsigned int n_phase_transitions;
260  };
261  }
262 
263  }
264 }
265 
266 #endif
void create_additional_named_outputs(MaterialModel::MaterialModelOutputs< dim > &out) const
double compute_partitioning_fraction(const double temperature) const
static void declare_parameters(ParameterHandler &prm)
void initialize_phase_function(std::shared_ptr< MaterialUtilities::PhaseFunction< dim >> &phase_function)
std::shared_ptr< MaterialUtilities::PhaseFunction< dim > > phase_function
Definition: compat.h:59
void calculate_reaction_terms(const typename Interface< dim >::MaterialModelInputs &in, const std::vector< double > &adiabatic_pressure, const std::vector< unsigned int > &phase_indices, const std::function< double(const double, const double, const double, const SymmetricTensor< 2, dim > &, const unsigned int, const double, const double)> &dislocation_viscosity, const std::function< double(const double, const double, const double, const double, const double, const unsigned int)> &diffusion_viscosity, const double min_eta, const double max_eta, typename Interface< dim >::MaterialModelOutputs &out) const
void fill_additional_outputs(const typename MaterialModel::MaterialModelInputs< dim > &in, const typename MaterialModel::MaterialModelOutputs< dim > &out, const std::vector< unsigned int > &phase_indices, const std::vector< double > &dislocation_viscosities, std::vector< std::unique_ptr< MaterialModel::AdditionalMaterialOutputs< dim >>> &additional_outputs) const