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  {
45  using namespace dealii;
46 
52  template <int dim>
54  {
55  public:
56  PlasticAdditionalOutputs(const unsigned int n_points);
57 
58  std::vector<double> get_nth_output(const unsigned int idx) const override;
59 
65  std::vector<double> cohesions;
66 
72  std::vector<double> friction_angles;
73 
77  std::vector<double> yield_stresses;
78 
83  std::vector<double> yielding;
84 
85  };
86 
91  {
95  std::vector<double> composition_viscosities;
96 
100  std::vector<bool> composition_yielding;
101 
105  std::vector<double> current_friction_angles;
106 
110  std::vector<double> current_cohesions;
111  };
112 
113  namespace Rheology
114  {
115 
116  template <int dim>
118  {
119  public:
123  ViscoPlastic();
124 
134  calculate_isostrain_viscosities ( const MaterialModel::MaterialModelInputs<dim> &in,
135  const unsigned int i,
136  const std::vector<double> &volume_fractions,
137  const std::vector<double> &phase_function_values = std::vector<double>(),
138  const std::vector<unsigned int> &n_phase_transitions_per_composition =
139  std::vector<unsigned int>()) const;
140 
150  void compute_viscosity_derivatives(const unsigned int point_index,
151  const std::vector<double> &volume_fractions,
152  const std::vector<double> &composition_viscosities,
155  const std::vector<double> &phase_function_values = std::vector<double>(),
156  const std::vector<unsigned int> &n_phase_transitions_per_composition =
157  std::vector<unsigned int>()) const;
158 
165  ComponentMask get_volumetric_composition_mask() const;
166 
170  static
171  void
172  declare_parameters (ParameterHandler &prm);
173 
181  void
182  parse_parameters (ParameterHandler &prm,
183  const std::unique_ptr<std::vector<unsigned int>> &expected_n_phases_per_composition = nullptr);
184 
189  void
190  create_plastic_outputs (MaterialModel::MaterialModelOutputs<dim> &out) const;
191 
197  void fill_plastic_outputs(const unsigned int point_index,
198  const std::vector<double> &volume_fractions,
199  const bool plastic_yielding,
202  const IsostrainViscosities &isostrain_viscosities) const;
203 
208 
213 
218 
223 
228 
229 
230  private:
231 
237 
243  std::vector<double> minimum_viscosity;
244  std::vector<double> maximum_viscosity;
245 
251  {
255  composite
256  } viscous_flow_law;
257 
263  {
265  drucker_prager
266  } yield_mechanism;
267 
274 
284 
293 
298  std::vector<double> exponents_stress_limiter;
299 
304 
310  std::unique_ptr<Rheology::FrankKamenetskii<dim>> frank_kamenetskii_rheology;
311 
316 
320  std::unique_ptr<Rheology::PeierlsCreep<dim>> peierls_creep;
321 
328 
333 
334  /*
335  * Object for computing plastic stresses, viscosities, and additional outputs
336  */
338 
339  /*
340  * Input parameters for the drucker prager plasticity.
341  */
343 
344  };
345  }
346  }
347 }
348 #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
Definition: compat.h:59
Rheology::CompositionalViscosityPrefactors< dim > compositional_viscosity_prefactors
Rheology::StrainDependent< dim > strain_rheology
Rheology::DiffusionCreep< dim > diffusion_creep
Definition: compat.h:42
std::unique_ptr< Rheology::PeierlsCreep< dim > > peierls_creep
Rheology::DruckerPragerParameters drucker_prager_parameters
Rheology::Elasticity< dim > elastic_rheology