Smith  0.1
Smith is an implicit thermal structural mechanics simulation code.
parameterized_solid_material.hpp
Go to the documentation of this file.
1 // Copyright (c) Lawrence Livermore National Security, LLC and
2 // other Smith Project Developers. See the top-level LICENSE file for
3 // details.
4 //
5 // SPDX-License-Identifier: (BSD-3-Clause)
6 
13 #pragma once
14 
17 
20 
26  using State = Empty;
27 
40  template <int dim, typename DispGradType, typename BulkType, typename ShearType>
42  const BulkType& DeltaK, const ShearType& DeltaG) const
43  {
44  constexpr auto I = Identity<dim>();
45  auto K = K0 + get<0>(DeltaK);
46  auto G = G0 + get<0>(DeltaG);
47  auto lambda = K - (2.0 / dim) * G;
48  auto epsilon = 0.5 * (transpose(du_dX) + du_dX);
49  return lambda * tr(epsilon) * I + 2.0 * G * epsilon;
50  }
51 
57  static constexpr int numParameters() { return 2; }
58 
59  double density;
60  double K0;
61  double G0;
62 };
63 
69  using State = Empty;
70 
83  template <int dim, typename DispGradType, typename BulkType, typename ShearType>
85  const BulkType& DeltaK, const ShearType& DeltaG) const
86  {
87  using std::log1p;
88  constexpr auto I = Identity<dim>();
89  auto K = K0 + get<0>(DeltaK);
90  auto G = G0 + get<0>(DeltaG);
91  auto lambda = K - (2.0 / dim) * G;
92  auto B_minus_I = du_dX * transpose(du_dX) + transpose(du_dX) + du_dX;
93  auto logJ = log1p(detApIm1(du_dX));
94 
95  // Kirchoff stress, in form that avoids cancellation error when F is near I
96  auto TK = lambda * logJ * I + G * B_minus_I;
97 
98  // Pull back to Piola
99  auto F = du_dX + I;
100  return dot(TK, inv(transpose(F)));
101  }
102 
108  static constexpr int numParameters() { return 2; }
109 
110  double density;
111  double K0;
112  double G0;
113 };
114 
120 template <typename... T>
122  using type = decltype((T{} + ...));
123 };
124 
127  static constexpr int dim = 3;
128  static constexpr double tol = 1e-10;
129 
130  double E;
131  double nu;
132  double density;
133 
135  struct State {
138  };
139 
141  template <typename DisplacementGradient, typename YieldStrength, typename SaturationStrength, typename StrainConstant>
142  auto operator()(State& state, const DisplacementGradient du_dX, const YieldStrength sigma_y,
143  const SaturationStrength sigma_sat, const StrainConstant strain_constant) const
144  {
145  // The output stress tensor should use dual numbers if any of the parameters are dual.
146  // This slightly ugly trick to accomplishes that by picking up any dual number types
147  // from the parameters into the dummy variable "one".
148  // Another possiblity would be to cast the results to the correct type at the end, which
149  // would avoid doing any unneccessary arithmetic with dual numbers.
151  T one;
152  one = 1.0;
153 
154  using std::sqrt;
155  constexpr auto I = Identity<dim>();
156  const double K = E / (3.0 * (1.0 - 2.0 * nu));
157  const double G = 0.5 * E / (1.0 + nu);
158  const double rel_tol = tol * get_value(sigma_y);
159 
160  // (i) elastic predictor
161  auto el_strain = sym(du_dX) - state.plastic_strain;
162  auto p = K * tr(el_strain);
163  auto s = 2.0 * G * dev(el_strain) * one; // multiply by "one" to get type correct for parameter derivatives
164  auto q = sqrt(1.5) * norm(s);
165 
166  // (ii) admissibility
167  const double eqps_old = state.accumulated_plastic_strain;
168  auto residual = [eqps_old, G, *this](auto delta_eqps, auto trial_mises, auto y0, auto ysat, auto e0) {
169  auto Y = this->flow_strength(eqps_old + delta_eqps, y0, ysat, e0);
170  return trial_mises - 3.0 * G * delta_eqps - Y;
171  };
172  if (residual(0.0, get_value(q), get_value(sigma_y), get_value(sigma_sat), get_value(strain_constant)) > rel_tol) {
173  // (iii) return mapping
174 
175  // Note the tolerance for convergence is the same as the tolerance for entering the return map.
176  // This ensures that if the constitutive update is called again with the updated internal
177  // variables, the return map won't be repeated.
178  ScalarSolverOptions opts{.xtol = 0, .rtol = rel_tol, .max_iter = 25};
179  double lower_bound = 0.0;
180  double upper_bound = (get_value(q) - flow_strength(eqps_old, get_value(sigma_y), get_value(sigma_sat),
181  get_value(strain_constant))) /
182  (3.0 * G);
183  auto [delta_eqps, status] =
184  solve_scalar_equation(residual, 0.0, lower_bound, upper_bound, opts, q, sigma_y, sigma_sat, strain_constant);
185 
186  auto Np = 1.5 * s / q;
187 
188  s = s - 2.0 * G * delta_eqps * Np;
189  state.accumulated_plastic_strain += get_value(delta_eqps);
190  state.plastic_strain += get_value(delta_eqps) * get_value(Np);
191  }
192 
193  return s + p * I;
194  }
195 
197  template <typename PlasticStrain, typename YieldStrength, typename SaturationStrength, typename StrainConstant>
198  auto flow_strength(const PlasticStrain accumulated_plastic_strain, const YieldStrength sigma_y,
199  const SaturationStrength sigma_sat, const StrainConstant strain_constant) const
200  {
201  using std::exp;
202  return sigma_sat - (sigma_sat - sigma_y) * exp(-accumulated_plastic_strain / strain_constant);
203  };
204 };
205 
206 } // namespace smith::solid_mechanics
#define SMITH_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc and does nothing on a host compi...
Definition: accelerator.hpp:38
Implementation of the quadrature-function-based functional enabling rapid development of FEM formulat...
SolidMechanics helper data types.
SMITH_HOST_DEVICE auto exp(dual< gradient_type > a)
implementation of exponential function for dual numbers
Definition: dual.hpp:381
constexpr SMITH_HOST_DEVICE auto detApIm1(const tensor< T, 2, 2 > &A)
computes det(A + I) - 1, where precision is not lost when the entries A_{ij} << 1
Definition: tensor.hpp:1319
auto solve_scalar_equation(const function &f, double x0, double lower_bound, double upper_bound, ScalarSolverOptions options, ParamTypes... params)
Solves a nonlinear scalar-valued equation and gives derivatives of solution to parameters.
constexpr SMITH_HOST_DEVICE auto inv(const isotropic_tensor< T, m, m > &I)
return the inverse of an isotropic tensor
SMITH_HOST_DEVICE auto sqrt(dual< gradient_type > x)
implementation of square root for dual numbers
Definition: dual.hpp:308
constexpr SMITH_HOST_DEVICE auto norm(const isotropic_tensor< T, m, m > &I)
compute the Frobenius norm (sqrt(tr(dot(transpose(I), I)))) of an isotropic tensor
constexpr SMITH_HOST_DEVICE auto transpose(const isotropic_tensor< T, m, m > &I)
return the transpose of an isotropic tensor
constexpr SMITH_HOST_DEVICE auto dev(const tensor< T, n, n > &A)
Calculates the deviator of a matrix (rank-2 tensor)
Definition: tensor.hpp:1193
constexpr SMITH_HOST_DEVICE auto get_value(const T &arg)
return the "value" part from a given type. For non-dual types, this is just the identity function
Definition: dual.hpp:445
constexpr SMITH_HOST_DEVICE auto sym(const isotropic_tensor< T, m, m > &I)
return the symmetric part of an isotropic tensor
constexpr SMITH_HOST_DEVICE auto tr(const isotropic_tensor< T, m, m > &I)
calculate the trace of an isotropic tensor
SMITH_HOST_DEVICE auto log1p(dual< gradient_type > a)
implementation of the natural logarithm of one plus the argument function for dual numbers
Definition: dual.hpp:397
constexpr SMITH_HOST_DEVICE auto dot(const isotropic_tensor< S, m, m > &I, const tensor< T, m, n... > &A)
dot product between an isotropic and (nonisotropic) tensor
The material and load types for the solid functional physics module.
see Nothing for a complete description of this class and when to use it
Settings for solve_scalar_equation.
double xtol
absolute tolerance on Newton correction
variables required to characterize the hysteresis response
J2 material with Voce hardening, with hardening parameters exposed as differentiable parameters.
static constexpr double tol
relative tolerance on residual for accepting return map solution
auto operator()(State &state, const DisplacementGradient du_dX, const YieldStrength sigma_y, const SaturationStrength sigma_sat, const StrainConstant strain_constant) const
calculate the first Piola stress, given the displacement gradient and previous material state
auto flow_strength(const PlasticStrain accumulated_plastic_strain, const YieldStrength sigma_y, const SaturationStrength sigma_sat, const StrainConstant strain_constant) const
Computes flow strength from Voce's hardening law.
static constexpr int numParameters()
The number of parameters in the model.
SMITH_HOST_DEVICE auto operator()(State &, const smith::tensor< DispGradType, dim, dim > &du_dX, const BulkType &DeltaK, const ShearType &DeltaG) const
stress calculation for a linear isotropic material model
SMITH_HOST_DEVICE auto operator()(State &, const smith::tensor< DispGradType, dim, dim > &du_dX, const BulkType &DeltaK, const ShearType &DeltaG) const
stress calculation for a NeoHookean material model
static constexpr int numParameters()
The number of parameters in the model.
Infers type resulting from algebraic expressions of a group of variables.
decltype((T{}+...)) type
type of the sum of the parameters
Arbitrary-rank tensor class.
Definition: tensor.hpp:28