ASPECT
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_visco_plastic_h
22 #define _aspect_material_model_rheology_visco_plastic_h
23 
24 #include <aspect/global.h>
38 
39 #include<deal.II/fe/component_mask.h>
40 
41 namespace aspect
42 {
43  namespace MaterialModel
44  {
50  template <int dim>
52  {
53  public:
54  PlasticAdditionalOutputs(const unsigned int n_points);
55 
56  std::vector<double> get_nth_output(const unsigned int idx) const override;
57 
63  std::vector<double> cohesions;
64 
70  std::vector<double> friction_angles;
71 
75  std::vector<double> yield_stresses;
76 
81  std::vector<double> yielding;
82 
83  };
84 
89  {
93  std::vector<double> composition_viscosities;
94 
98  std::vector<bool> composition_yielding;
99 
103  std::vector<double> current_friction_angles;
104 
108  std::vector<double> current_cohesions;
109  };
110 
111  namespace Rheology
112  {
113 
114  template <int dim>
116  {
117  public:
121  ViscoPlastic();
122 
132  calculate_isostrain_viscosities ( const MaterialModel::MaterialModelInputs<dim> &in,
133  const unsigned int i,
134  const std::vector<double> &volume_fractions,
135  const std::vector<double> &phase_function_values = std::vector<double>(),
136  const std::vector<unsigned int> &n_phase_transitions_per_composition =
137  std::vector<unsigned int>()) const;
138 
148  void compute_viscosity_derivatives(const unsigned int point_index,
149  const std::vector<double> &volume_fractions,
150  const std::vector<double> &composition_viscosities,
153  const std::vector<double> &phase_function_values = std::vector<double>(),
154  const std::vector<unsigned int> &n_phase_transitions_per_composition =
155  std::vector<unsigned int>()) const;
156 
163  ComponentMask get_volumetric_composition_mask() const;
164 
168  static
169  void
170  declare_parameters (ParameterHandler &prm);
171 
179  void
180  parse_parameters (ParameterHandler &prm,
181  const std::unique_ptr<std::vector<unsigned int>> &expected_n_phases_per_composition = nullptr);
182 
187  void
188  create_plastic_outputs (MaterialModel::MaterialModelOutputs<dim> &out) const;
189 
195  void fill_plastic_outputs(const unsigned int point_index,
196  const std::vector<double> &volume_fractions,
197  const bool plastic_yielding,
200  const IsostrainViscosities &isostrain_viscosities) const;
201 
206 
211 
216 
221 
226 
227 
228  private:
229 
235 
241  std::vector<double> minimum_viscosity;
242  std::vector<double> maximum_viscosity;
243 
249  {
253  composite
254  } viscous_flow_law;
255 
261  {
263  drucker_prager
264  } yield_mechanism;
265 
272 
282 
291 
296  std::vector<double> exponents_stress_limiter;
297 
302 
308  std::unique_ptr<Rheology::FrankKamenetskii<dim>> frank_kamenetskii_rheology;
309 
314 
318  std::unique_ptr<Rheology::PeierlsCreep<dim>> peierls_creep;
319 
326 
331 
332  /*
333  * Object for computing plastic stresses, viscosities, and additional outputs
334  */
336 
337  /*
338  * Input parameters for the drucker prager plasticity.
339  */
341 
342  };
343  }
344  }
345 }
346 #endif
MaterialUtilities::CompositionalAveragingOperation viscosity_averaging
Rheology::FrictionModels< dim > friction_models
std::unique_ptr< Rheology::FrankKamenetskii< dim > > frank_kamenetskii_rheology
Rheology::DruckerPrager< dim > drucker_prager_plasticity
void declare_parameters(ParameterHandler &prm)
Rheology::ConstantViscosityPrefactors< dim > constant_viscosity_prefactors
Rheology::DislocationCreep< dim > dislocation_creep
PlasticAdditionalOutputs(const unsigned int n_points)
Definition: compat.h:59
Rheology::CompositionalViscosityPrefactors< dim > compositional_viscosity_prefactors
Rheology::StrainDependent< dim > strain_rheology
Rheology::DiffusionCreep< dim > diffusion_creep
std::unique_ptr< Rheology::PeierlsCreep< dim > > peierls_creep
Rheology::DruckerPragerParameters drucker_prager_parameters
Rheology::Elasticity< dim > elastic_rheology
std::vector< double > get_nth_output(const unsigned int idx) const override