38 template <
typename T,
int dim>
41 auto I = Identity<dim>();
42 auto lambda =
K - (2.0 / 3.0) *
G;
43 auto epsilon = 0.5 * (
transpose(du_dX) + du_dX);
44 return lambda *
tr(epsilon) * I + 2.0 *
G * epsilon;
55 template <
typename T,
int dim>
74 template <
typename T,
int dim>
77 static constexpr
auto I = Identity<dim>();
82 const auto S =
K *
tr(E) * I + 2.0 *
G *
dev(E);
109 template <
typename T,
int dim>
113 constexpr
auto I = Identity<dim>();
114 auto lambda =
K - (2.0 / 3.0) *
G;
119 auto TK = lambda * logJ * I +
G * B_minus_I;
142 template <
typename T,
int dim>
146 constexpr
auto I = Identity<dim>();
149 auto Jm13 =
pow(J, -1.0 / 3.0);
150 auto F_bar = Jm13 * F;
151 auto Pdev =
G * Jm13 * (F_bar - 1.0 / 3.0 *
inner(F_bar, F_bar) *
inv(
transpose(F_bar)));
169 template <
typename T>
170 auto overstress(
double viscosity, T accumulated_plastic_strain_rate)
172 return viscosity * accumulated_plastic_strain_rate;
192 template <
typename T1,
typename T2>
193 auto operator()(T1 accumulated_plastic_strain, T2 accumulated_plastic_strain_rate)
const
217 template <
typename T1,
typename T2>
218 auto operator()(T1 accumulated_plastic_strain, T2 accumulated_plastic_strain_rate)
const
221 return sigma_y *
pow(1.0 + accumulated_plastic_strain /
eps0, 1.0 /
n) +
246 template <
typename T1,
typename T2>
247 auto operator()(T1 accumulated_plastic_strain, T2 accumulated_plastic_strain_rate)
const
250 auto& epsilon_p = accumulated_plastic_strain;
252 auto& epsilon_p_dot = accumulated_plastic_strain_rate;
259 template <
typename HardeningType>
261 static constexpr
int dim = 3;
262 static constexpr
double tol = 1e-10;
277 template <
typename T>
281 constexpr
auto I = Identity<dim>();
282 const double K =
E / (3.0 * (1.0 - 2.0 *
nu));
283 const double G = 0.5 *
E / (1.0 +
nu);
287 auto p = K *
tr(el_strain);
288 auto s = 2.0 * G *
dev(el_strain);
290 auto eta = s - sigma_b;
295 auto residual = [eqps_old, G, dt, *
this](
auto delta_eqps,
auto trial_q) {
296 auto eqps_dot = dt > 0.0 ? delta_eqps / dt : 0.0 * delta_eqps;
297 return trial_q - (3.0 * G +
Hk) * delta_eqps - this->
hardening(eqps_old + delta_eqps, eqps_dot);
306 double lower_bound = 0.0;
310 auto Np = 1.5 * eta / q;
312 s = s - 2.0 * G * delta_eqps * Np;
322 template <
typename HardeningType>
324 static constexpr
int dim = 3;
325 static constexpr
double tol = 1e-10;
339 template <
typename T>
343 constexpr
auto I = Identity<dim>();
344 const double K =
E / (3.0 * (1.0 - 2.0 *
nu));
345 const double G = 0.5 *
E / (1.0 +
nu);
354 auto s = 2.0 * G *
dev(Ee);
359 auto residual = [eqps_old, G, dt, *
this](
auto delta_eqps,
auto trial_mises) {
360 auto eqps_dot = dt > 0.0 ? delta_eqps / dt : 0.0 * delta_eqps;
361 return trial_mises - 3.0 * G * delta_eqps - this->
hardening(eqps_old + delta_eqps, eqps_dot);
370 double lower_bound = 0.0;
371 double upper_bound = (q_val -
hardening(eqps_old, 0.0)) / (3.0 * G);
376 if (!status.converged) {
377 SLIC_WARNING(
"J2 solve failed to converge.");
380 auto Np = 1.5 * s / q;
382 s = s - 2.0 * G * delta_eqps * Np;
383 auto A =
exp_symm(-delta_eqps * Np);
418 template <
typename T1,
typename T2,
int dim>
421 return transpose(
dot(
inv(displacement_gradient + Identity<dim>()), kirchhoff_stress));
435 template <
typename T1,
typename T2,
int dim>
438 auto kirchhoff_stress =
det(displacement_gradient + Identity<dim>()) * cauchy_stress;
455 template <
typename T>
474 template <
typename T>
502 template <
typename T,
int dim,
typename Density>
506 constexpr
auto I = Identity<dim>();
507 auto lambda =
K - (2.0 / 3.0) *
G;
512 auto TK = lambda * logJ * I +
G * B_minus_I;
520 template <
typename Density>
#define SMITH_HOST_DEVICE
Macro that evaluates to __host__ __device__ when compiling with nvcc or amdclang and does nothing on ...
Implementation of the quadrature-function-based functional enabling rapid development of FEM formulat...
SolidMechanics helper data types.
auto overstress(double viscosity, T accumulated_plastic_strain_rate)
auto CauchyToPiola(const tensor< T1, dim, dim > &cauchy_stress, const tensor< T2, dim, dim > &displacement_gradient)
Transform the Cauchy stress to the Piola stress.
auto greenStrain(const tensor< T, dim, dim > &grad_u)
Compute Green's strain from the displacement gradient.
auto KirchhoffToPiola(const tensor< T1, dim, dim > &kirchhoff_stress, const tensor< T2, dim, dim > &displacement_gradient)
Transform the Kirchhoff stress to the Piola stress.
SMITH_HOST_DEVICE auto exp(dual< gradient_type > a)
implementation of exponential function for dual numbers
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
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
auto log_symm(tensor< T, 3, 3 > A)
Logarithm of a symmetric matrix.
SMITH_HOST_DEVICE auto sqrt(dual< gradient_type > x)
implementation of square root for dual numbers
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)
SMITH_HOST_DEVICE auto pow(dual< gradient_type > a, dual< gradient_type > b)
implementation of a (dual) raised to the b (dual) power
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
constexpr SMITH_HOST_DEVICE auto sym(const isotropic_tensor< T, m, m > &I)
return the symmetric part of an isotropic tensor
auto exp_symm(tensor< T, 3, 3 > A)
Exponential of a symmetric matrix.
constexpr SMITH_HOST_DEVICE auto det(const isotropic_tensor< T, m, m > &I)
compute the determinant of an isotropic tensor
constexpr SMITH_HOST_DEVICE auto inner(const dual< S > &A, const dual< T > &B)
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
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
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
Constant body force model.
tensor< double, dim > force_
The constant body force.
SMITH_HOST_DEVICE tensor< double, dim > operator()(const tensor< T, dim > &, const double) const
Evaluation function for the constant body force model.
Constant traction boundary condition model.
SMITH_HOST_DEVICE tensor< double, dim > operator()(const tensor< T, dim > &, const tensor< T, dim > &, const double) const
Evaluation function for the constant traction model.
tensor< double, dim > traction_
The constant traction.
variables required to characterize the hysteresis response
double accumulated_plastic_strain
uniaxial equivalent plastic strain
tensor< double, dim, dim > plastic_strain
plastic strain
J2 material with nonlinear isotropic hardening and linear kinematic hardening.
static constexpr double tol
relative tolerance on residual mag to judge convergence of return map
auto operator()(State &state, double dt, const tensor< T, dim, dim > &du_dX) const
calculate the first Piola stress, given the displacement gradient and previous material state
double nu
Poisson's ratio.
double Hk
Kinematic hardening modulus.
double density
Mass density.
static constexpr int dim
spatial dimension
HardeningType hardening
Flow stress hardening model.
variables required to characterize the hysteresis response
tensor< double, dim, dim > Fpinv
inverse of plastic distortion tensor
double accumulated_plastic_strain
uniaxial equivalent plastic strain
Finite deformation version of J2 material with nonlinear isotropic hardening.
double density
mass density
static constexpr int dim
spatial dimension
HardeningType hardening
Flow stress hardening model.
static constexpr double tol
relative tolerance on residual mag to judge convergence of return map
auto operator()(State &state, double dt, const tensor< T, dim, dim > &du_dX) const
calculate the first Piola stress, given the displacement gradient and previous material state
double nu
Poisson's ratio.
Linear isotropic hardening law.
double Hi
Isotropic hardening modulus.
double eta
viscosity for linear rate sensitivity
auto operator()(T1 accumulated_plastic_strain, T2 accumulated_plastic_strain_rate) const
Computes the flow stress.
double sigma_y
yield strength
Linear isotropic elasticity material model.
double density
mass density
SMITH_HOST_DEVICE auto operator()(State &, const tensor< T, dim, dim > &du_dX) const
stress calculation for a linear isotropic material model
Neo-Hookean material version with additive split of deviatoric and volumetric behavior.
SMITH_HOST_DEVICE auto operator()(State &, const tensor< T, dim, dim > &du_dX) const
Piola stress calculation.
double density
mass density
Neo-Hookean material model This struct differs in style relative to the older materials as it needs t...
SMITH_HOST_DEVICE auto density(const Density &density) const
interpolates density field
SMITH_HOST_DEVICE auto pkStress(State &, const tensor< T, dim, dim > &du_dX, const Density &) const
stress calculation for a NeoHookean material model
Neo-Hookean material model.
SMITH_HOST_DEVICE auto operator()(State &, const tensor< T, dim, dim > &du_dX) const
stress calculation for a NeoHookean material model
double density
mass density
Power-law isotropic hardening law.
double eta
viscosity for linear rate sensitivity
double sigma_y
yield strength
double eps0
reference value of accumulated plastic strain
double n
hardening index in reciprocal form
auto operator()(T1 accumulated_plastic_strain, T2 accumulated_plastic_strain_rate) const
Computes the flow stress.
St. Venant Kirchhoff hyperelastic model.
auto operator()(State &, const tensor< T, dim, dim > &grad_u) const
stress calculation for a St. Venant Kirchhoff material model
Voce's isotropic hardening law.
double strain_constant
The constant dictating how fast the exponential decays.
double sigma_y
yield strength
double sigma_sat
saturation value of flow strength
auto operator()(T1 accumulated_plastic_strain, T2 accumulated_plastic_strain_rate) const
Computes the flow stress.
double eta
viscosity for linear rate sensitivity
Arbitrary-rank tensor class.