https://mooseframework.inl.gov
LagrangianStressDivergenceBaseS.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 
14 {
16 
17  params.addRequiredParam<unsigned int>("component", "Which direction this kernel acts in");
18  params.addRequiredCoupledVar("displacements", "The displacement components");
19 
20  params.addParam<bool>("large_kinematics", false, "Use large displacement kinematics");
21  params.addParam<bool>("stabilize_strain", false, "Average the volumetric strains");
22 
23  params.addParam<std::string>("base_name", "Material property base name");
24 
25  params.addCoupledVar("temperature",
26  "The name of the temperature variable used in the "
27  "ComputeThermalExpansionEigenstrain. (Not required for "
28  "simulations without temperature coupling.)");
29 
30  params.addParam<std::vector<MaterialPropertyName>>(
31  "eigenstrain_names",
32  {},
33  "List of eigenstrains used in the strain calculation. Used for computing their derivatives "
34  "for off-diagonal Jacobian terms.");
35 
36  params.addCoupledVar("out_of_plane_strain",
37  "The out-of-plane strain variable for weak plane stress formulation.");
38 
39  return params;
40 }
41 
44  _large_kinematics(getParam<bool>("large_kinematics")),
45  _stabilize_strain(getParam<bool>("stabilize_strain")),
46  _base_name(isParamValid("base_name") ? getParam<std::string>("base_name") + "_" : ""),
47  _alpha(getParam<unsigned int>("component")),
48  _ndisp(coupledComponents("displacements")),
49  _disp_nums(_ndisp),
50  _avg_grad_trial(_ndisp),
51  _F_ust(
52  getMaterialPropertyByName<RankTwoTensor>(_base_name + "unstabilized_deformation_gradient")),
53  _F_avg(getMaterialPropertyByName<RankTwoTensor>(_base_name + "average_deformation_gradient")),
54  _f_inv(getMaterialPropertyByName<RankTwoTensor>(_base_name +
55  "inverse_incremental_deformation_gradient")),
56  _F_inv(getMaterialPropertyByName<RankTwoTensor>(_base_name + "inverse_deformation_gradient")),
57  _F(getMaterialPropertyByName<RankTwoTensor>(_base_name + "deformation_gradient")),
58  _temperature(isCoupled("temperature") ? getVar("temperature", 0) : nullptr),
59  _out_of_plane_strain(isCoupled("out_of_plane_strain") ? getVar("out_of_plane_strain", 0)
60  : nullptr)
61 {
62  // Do the vector coupling of the displacements
63  for (unsigned int i = 0; i < _ndisp; i++)
64  _disp_nums[i] = coupled("displacements", i);
65 
66  // We need to use identical discretizations for all displacement components
67  auto order_x = getVar("displacements", 0)->order();
68  for (unsigned int i = 1; i < _ndisp; i++)
69  {
70  if (getVar("displacements", i)->order() != order_x)
71  mooseError("The Lagrangian StressDivergence kernels require equal "
72  "order interpolation for all displacements.");
73  }
74 
75  // fetch eigenstrain derivatives
76  const auto nvar = _coupled_moose_vars.size();
77  _deigenstrain_dargs.resize(nvar);
78  for (std::size_t i = 0; i < nvar; ++i)
79  for (auto eigenstrain_name : getParam<std::vector<MaterialPropertyName>>("eigenstrain_names"))
80  _deigenstrain_dargs[i].push_back(&getMaterialPropertyDerivative<RankTwoTensor>(
81  eigenstrain_name, _coupled_moose_vars[i]->name()));
82 }
83 
84 void
86 {
87  // Skip if we are not doing stabilization
88  if (!_stabilize_strain)
89  return;
90 
91  // We need the gradients of shape functions in the reference frame
92  _fe_problem.prepareShapes(_var.number(), _tid);
93  _avg_grad_trial[_alpha].resize(_phi.size());
95 }
96 
97 void
99 {
100  // Skip if we are not doing stabilization
101  if (!_stabilize_strain)
102  return;
103 
104  for (auto beta : make_range(_ndisp))
105  if (jvar == _disp_nums[beta])
106  {
107  // We need the gradients of shape functions in the reference frame
108  _fe_problem.prepareShapes(jvar, _tid);
109  _avg_grad_trial[beta].resize(_phi.size());
111  }
112 }
113 
114 Real
116 {
118 }
119 
120 Real
122 {
123  // Bail if jvar not coupled
124  if (getJvarMap()[jvar] < 0)
125  return 0.0;
126 
127  // Off diagonal terms for other displacements
128  for (auto beta : make_range(_ndisp))
129  if (jvar == _disp_nums[beta])
131 
132  // Off diagonal temperature term due to eigenstrain
133  if (_temperature && jvar == _temperature->number())
135 
136  // Off diagonal term due to weak plane stress
139 
140  return 0;
141 }
const MooseVariable * _temperature
Temperature, if provided. This is used only to get the trial functions.
std::vector< unsigned int > _disp_nums
The displacement numbers.
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void mooseError(Args &&... args)
unsigned int number() const
void addRequiredParam(const std::string &name, const std::string &doc_string)
const MooseVariable * _out_of_plane_strain
Out-of-plane strain, if provided.
const unsigned int _ndisp
Total number of displacements/size of residual vector.
const std::string name
Definition: Setup.h:20
virtual Real computeQpOffDiagJacobian(unsigned int jvar) override
void addCoupledVar(const std::string &name, const std::string &doc_string)
virtual Real computeQpJacobianOutOfPlaneStrain()=0
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
const bool _stabilize_strain
If true calculate the deformation gradient derivatives for F_bar.
std::vector< std::vector< const MaterialProperty< RankTwoTensor > * > > _deigenstrain_dargs
Eigenstrain derivatives wrt generate coupleds.
virtual Real computeQpJacobianDisplacement(unsigned int alpha, unsigned int beta)=0
LagrangianStressDivergenceBaseS(const InputParameters &parameters)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
virtual Real computeQpJacobianTemperature(unsigned int cvar)=0
static InputParameters validParams()
IntRange< T > make_range(T beg, T end)
virtual void precalculateJacobianDisplacement(unsigned int component)=0
Prepare the average shape function gradients for stabilization.
void ErrorVector unsigned int
virtual void precalculateOffDiagJacobian(unsigned int jvar) override