https://mooseframework.inl.gov
TotalLagrangianStressDivergenceBaseS.C
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
11 
12 template <class G>
15 {
17  // This kernel requires use_displaced_mesh to be off
18  params.suppressParameter<bool>("use_displaced_mesh");
19  return params;
20 }
21 
22 template <class G>
24  const InputParameters & parameters)
25  : LagrangianStressDivergenceBaseS(parameters),
26  _pk1(getMaterialPropertyByName<RankTwoTensor>(_base_name + "pk1_stress")),
27  _dpk1(getMaterialPropertyByName<RankFourTensor>(_base_name + "pk1_jacobian"))
28 {
29 }
30 
31 template <class G>
34 {
35  // F-bar doesn't modify the test function
36  return G::gradOp(component, _grad_test[_i][_qp], _test[_i][_qp], _q_point[_qp]);
37 }
38 
39 template <class G>
42 {
43  return _stabilize_strain ? gradTrialStabilized(component) : gradTrialUnstabilized(component);
44 }
45 
46 template <class G>
49 {
50  // Without F-bar stabilization, simply return the gradient of the trial functions
51  return G::gradOp(component, _grad_phi[_j][_qp], _phi[_j][_qp], _q_point[_qp]);
52 }
53 
54 template <class G>
57 {
58  // The base unstabilized trial function gradient
59  const auto Gb = G::gradOp(component, _grad_phi[_j][_qp], _phi[_j][_qp], _q_point[_qp]);
60  // The average trial function gradient
61  const auto Ga = _avg_grad_trial[component][_j];
62 
63  // The F-bar stabilization depends on kinematics
64  if (_large_kinematics)
65  {
66  // Horrible thing, see the documentation for how we get here
67  const Real dratio = std::pow(_F_avg[_qp].det() / _F_ust[_qp].det(), 1.0 / 3.0);
68  const Real fact = (_F_avg[_qp].inverse().transpose().doubleContraction(Ga) -
69  _F_ust[_qp].inverse().transpose().doubleContraction(Gb)) /
70  3.0;
71  return dratio * (Gb + fact * _F_ust[_qp]);
72  }
73 
74  // The small kinematics modification is linear
75  return Gb + (Ga.trace() - Gb.trace()) / 3.0 * RankTwoTensor::Identity();
76 }
77 
78 template <class G>
79 void
81 {
82  // For total Lagrangian, the averaging is taken on the reference frame regardless of geometric
83  // nonlinearity. Convenient!
84  for (auto j : make_range(_phi.size()))
86  [this, component, j](unsigned int qp)
87  { return G::gradOp(component, _grad_phi[j][qp], _phi[j][qp], _q_point[qp]); },
88  _JxW,
89  _coord);
90 }
91 
92 template <class G>
93 Real
95 {
96  return gradTest(_alpha).doubleContraction(_pk1[_qp]);
97 }
98 
99 template <class G>
100 Real
102  unsigned int beta)
103 {
104  // J_{alpha beta} = phi^alpha_{i, J} T_{iJkL} G^beta_{kL}
105  return gradTest(alpha).doubleContraction(_dpk1[_qp] * gradTrial(beta));
106 }
107 
108 template <class G>
109 Real
111 {
112  usingTensorIndices(i_, j_, k_, l_);
113  // Multiple eigenstrains may depend on the same coupled var
114  RankTwoTensor total_deigen;
115  for (const auto deigen_darg : _deigenstrain_dargs[cvar])
116  total_deigen += (*deigen_darg)[_qp];
117 
118  const auto A = _f_inv[_qp].inverse();
119  const auto B = _F_inv[_qp].inverse();
120  const auto U = 0.5 * (A.template times<i_, k_, l_, j_>(B) + A.template times<i_, l_, k_, j_>(B));
121 
122  return -(_dpk1[_qp] * U * total_deigen).doubleContraction(gradTest(_alpha)) *
123  _temperature->phi()[_j][_qp];
124 }
125 
126 template <class G>
127 Real
129 {
130  return _dpk1[_qp].contractionKl(2, 2, gradTest(_alpha)) * _out_of_plane_strain->phi()[_j][_qp];
131 }
132 
virtual void precalculateJacobianDisplacement(unsigned int component) override
Prepare the average shape function gradients for stabilization.
static const std::string component
Definition: NS.h:153
TotalLagrangianStressDivergenceBaseS(const InputParameters &parameters)
virtual Real computeQpJacobianDisplacement(unsigned int alpha, unsigned int beta) override
Enforce equilibrium with a total Lagrangian formulation.
static RankTwoTensorTempl Identity()
void suppressParameter(const std::string &name)
virtual RankTwoTensor gradTrialUnstabilized(unsigned int component)
The unstabilized trial function gradient.
virtual RankTwoTensor gradTrialStabilized(unsigned int component)
The stabilized trial function gradient.
virtual RankTwoTensor gradTest(unsigned int component) override
virtual Real computeQpJacobianTemperature(unsigned int cvar) override
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::string alpha
Definition: NS.h:134
IntRange< T > make_range(T beg, T end)
auto elementAverage(const Functor &f, const MooseArray< Real > &JxW, const MooseArray< Real > &coord)
Base class of the "Lagrangian" kernel system.
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
MooseUnits pow(const MooseUnits &, int)
virtual RankTwoTensor gradTrial(unsigned int component) override