https://mooseframework.inl.gov
Public Types | Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Attributes | List of all members
HomogenizedTotalLagrangianStressDivergence Class Reference

Total Lagrangian formulation with all homogenization terms (one disp_xyz field and macro_gradient scalar) More...

#include <HomogenizedTotalLagrangianStressDivergence.h>

Inheritance diagram for HomogenizedTotalLagrangianStressDivergence:
[legend]

Public Types

typedef std::vector< intJvarMap
 

Public Member Functions

 HomogenizedTotalLagrangianStressDivergence (const InputParameters &parameters)
 
virtual std::set< std::string > additionalROVariables () override
 Inform moose that this kernel covers the constraint scalar variable. More...
 
virtual void initialSetup () override
 
template<>
void initialSetup ()
 
template<>
void initialSetup ()
 
template<>
void initialSetup ()
 
virtual void computeOffDiagJacobian (unsigned int jvar) override
 
unsigned int mapJvarToCvar (unsigned int jvar)
 
int mapJvarToCvar (unsigned int jvar, const JvarMap &jvar_map)
 
bool mapJvarToCvar (unsigned int jvar, unsigned int &cvar)
 
const JvarMapgetJvarMap ()
 
const JvarMapgetParameterJvarMap (std::string parameter_name)
 

Static Public Member Functions

static InputParameters validParams ()
 
static InputParameters baseParams ()
 

Protected Member Functions

virtual void computeScalarResidual () override
 Method for computing the scalar part of residual for _kappa. More...
 
virtual void computeScalarJacobian () override
 Method for computing the scalar variable part of Jacobian for d-_kappa-residual / d-_kappa. More...
 
virtual void computeScalarOffDiagJacobian (const unsigned int jvar_num) override
 Method for computing an off-diagonal jacobian component d-_kappa-residual / d-jvar. More...
 
virtual Real computeScalarQpOffDiagJacobian (const unsigned int jvar_num) override
 Method for computing an off-diagonal jacobian component at quadrature points. More...
 
virtual void computeOffDiagJacobianScalarLocal (const unsigned int svar_num) override
 Method for computing an off-diagonal jacobian component d-_var-residual / d-svar. More...
 
virtual Real computeQpOffDiagJacobianScalar (const unsigned int svar_num) override
 Method for computing d-_var-residual / d-svar at quadrature points. More...
 
const Homogenization::ConstraintMapcmap () const
 Get the constraint map. More...
 
virtual RankTwoTensor gradTest (unsigned int component) override
 
virtual RankTwoTensor gradTrial (unsigned int component) override
 
virtual void precalculateJacobianDisplacement (unsigned int component) override
 Prepare the average shape function gradients for stabilization. More...
 
virtual Real computeQpResidual () override
 
virtual Real computeQpJacobianDisplacement (unsigned int alpha, unsigned int beta) override
 
virtual Real computeQpJacobianTemperature (unsigned int cvar) override
 
virtual Real computeQpJacobianOutOfPlaneStrain () override
 
virtual void precalculateJacobian () override
 
virtual void precalculateOffDiagJacobian (unsigned int jvar) override
 
virtual Real computeQpJacobian () override
 
virtual Real computeQpOffDiagJacobian (unsigned int jvar) override
 

Protected Attributes

const MaterialProperty< RankTwoTensor > & _pk1
 The 1st Piola-Kirchhoff stress. More...
 
const MaterialProperty< RankFourTensor > & _dpk1
 The derivative of the PK1 stress with respect to the deformation gradient. More...
 
const bool _large_kinematics
 If true use large deformation kinematics. More...
 
const bool _stabilize_strain
 If true calculate the deformation gradient derivatives for F_bar. More...
 
const std::string _base_name
 Prepend to the material properties. More...
 
const unsigned int _alpha
 Which component of the vector residual this kernel is responsible for. More...
 
const unsigned int _ndisp
 Total number of displacements/size of residual vector. More...
 
std::vector< unsigned int_disp_nums
 The displacement numbers. More...
 
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
 
const MaterialProperty< RankTwoTensor > & _F_ust
 The unmodified deformation gradient. More...
 
const MaterialProperty< RankTwoTensor > & _F_avg
 The element-average deformation gradient. More...
 
const MaterialProperty< RankTwoTensor > & _f_inv
 The inverse increment deformation gradient. More...
 
const MaterialProperty< RankTwoTensor > & _F_inv
 The inverse deformation gradient. More...
 
const MaterialProperty< RankTwoTensor > & _F
 The actual (stabilized) deformation gradient. More...
 
const MooseVariable_temperature
 Temperature, if provided. This is used only to get the trial functions. More...
 
const MooseVariable_out_of_plane_strain
 Out-of-plane strain, if provided. More...
 
std::vector< std::vector< const MaterialProperty< RankTwoTensor > * > > _deigenstrain_dargs
 Eigenstrain derivatives wrt generate coupleds. More...
 
const unsigned int _n_args
 

Private Attributes

unsigned int _m
 Indices for off-diagonal Jacobian components. More...
 
unsigned int _n
 
Homogenization::ConstraintType _ctype = Homogenization::ConstraintType::None
 Type of current homogenization constraint. More...
 

Detailed Description

Total Lagrangian formulation with all homogenization terms (one disp_xyz field and macro_gradient scalar)

Definition at line 17 of file HomogenizedTotalLagrangianStressDivergence.h.

Constructor & Destructor Documentation

◆ HomogenizedTotalLagrangianStressDivergence()

HomogenizedTotalLagrangianStressDivergence::HomogenizedTotalLagrangianStressDivergence ( const InputParameters parameters)

Member Function Documentation

◆ additionalROVariables()

std::set< std::string > HomogenizedTotalLagrangianStressDivergence::additionalROVariables ( )
overridevirtual

Inform moose that this kernel covers the constraint scalar variable.

Definition at line 37 of file HomogenizedTotalLagrangianStressDivergence.C.

38 {
39  // Add the scalar variable to the list of variables that this kernel contributes to
40  std::set<std::string> vars = TotalLagrangianStressDivergence::additionalROVariables();
41  vars.insert(_kappa_var_ptr->name());
42  return vars;
43 }
char ** vars

◆ baseParams()

template<class G >
InputParameters TotalLagrangianStressDivergenceBase< G >::baseParams ( )
staticinherited

Definition at line 14 of file TotalLagrangianStressDivergenceBase.C.

Referenced by TotalLagrangianStressDivergenceBase< G >::validParams().

15 {
17  // This kernel requires use_displaced_mesh to be off
18  params.suppressParameter<bool>("use_displaced_mesh");
19  return params;
20 }
void suppressParameter(const std::string &name)

◆ cmap()

Get the constraint map.

Definition at line 46 of file HomogenizationInterface.h.

Referenced by computeOffDiagJacobianScalarLocal(), computeScalarJacobian(), computeScalarOffDiagJacobian(), and computeScalarResidual().

46 { return _cmap; }
Homogenization::ConstraintMap _cmap
Constraint map.

◆ computeOffDiagJacobianScalarLocal()

void HomogenizedTotalLagrangianStressDivergence::computeOffDiagJacobianScalarLocal ( const unsigned int  svar_num)
overrideprotectedvirtual

Method for computing an off-diagonal jacobian component d-_var-residual / d-svar.

svar is looped over all scalar variables, which herein is just _kappa

Definition at line 187 of file HomogenizedTotalLagrangianStressDivergence.C.

189 {
190  // Just in case, skip any other scalar variables
191  if (svar_num != _kappa_var)
192  return;
193 
194  _local_ke.resize(_test.size(), _k_order);
195 
196  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
197  {
198  unsigned int l = 0;
199  const auto dV = _JxW[_qp] * _coord[_qp];
200  for (const auto & [indices, constraint] : cmap())
201  {
202  // copy constraint indices to protected variables to pass to Qp routine
203  std::tie(_m, _n) = indices;
204  _ctype = constraint.first;
205  initScalarQpJacobian(svar_num);
206  for (_i = 0; _i < _test.size(); _i++)
207  _local_ke(_i, l) += dV * computeQpOffDiagJacobianScalar(svar_num);
208  l++;
209  }
210  }
211 
212  addJacobian(
213  _assembly, _local_ke, _var.dofIndices(), _kappa_var_ptr->dofIndices(), _var.scalingFactor());
214 }
virtual Real computeQpOffDiagJacobianScalar(const unsigned int svar_num) override
Method for computing d-_var-residual / d-svar at quadrature points.
unsigned int _m
Indices for off-diagonal Jacobian components.
Homogenization::ConstraintType _ctype
Type of current homogenization constraint.
const Homogenization::ConstraintMap & cmap() const
Get the constraint map.

◆ computeQpJacobian()

Real LagrangianStressDivergenceBase::computeQpJacobian ( )
overrideprotectedvirtualinherited

Reimplemented in TotalLagrangianWeakPlaneStress.

Definition at line 115 of file LagrangianStressDivergenceBase.C.

116 {
118 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
virtual Real computeQpJacobianDisplacement(unsigned int alpha, unsigned int beta)=0

◆ computeQpJacobianDisplacement()

template<class G >
Real TotalLagrangianStressDivergenceBase< G >::computeQpJacobianDisplacement ( unsigned int  alpha,
unsigned int  beta 
)
overrideprotectedvirtualinherited

Implements LagrangianStressDivergenceBase.

Definition at line 101 of file TotalLagrangianStressDivergenceBase.C.

103 {
104  // J_{alpha beta} = phi^alpha_{i, J} T_{iJkL} G^beta_{kL}
105  return gradTest(alpha).doubleContraction(_dpk1[_qp] * gradTrial(beta));
106 }
virtual RankTwoTensor gradTest(unsigned int component) override
Real doubleContraction(const RankTwoTensorTempl< Real > &a) const
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.
static const std::string alpha
Definition: NS.h:134
virtual RankTwoTensor gradTrial(unsigned int component) override

◆ computeQpJacobianOutOfPlaneStrain()

template<class G >
Real TotalLagrangianStressDivergenceBase< G >::computeQpJacobianOutOfPlaneStrain ( )
overrideprotectedvirtualinherited

Implements LagrangianStressDivergenceBase.

Definition at line 128 of file TotalLagrangianStressDivergenceBase.C.

129 {
130  return _dpk1[_qp].contractionKl(2, 2, gradTest(_alpha)) * _out_of_plane_strain->phi()[_j][_qp];
131 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
const MooseVariable * _out_of_plane_strain
Out-of-plane strain, if provided.
virtual RankTwoTensor gradTest(unsigned int component) override
const FieldVariablePhiValue & phi() const override
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.

◆ computeQpJacobianTemperature()

template<class G >
Real TotalLagrangianStressDivergenceBase< G >::computeQpJacobianTemperature ( unsigned int  cvar)
overrideprotectedvirtualinherited

Implements LagrangianStressDivergenceBase.

Definition at line 110 of file TotalLagrangianStressDivergenceBase.C.

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 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
const MaterialProperty< RankTwoTensor > & _F_inv
The inverse deformation gradient.
virtual RankTwoTensor gradTest(unsigned int component) override
const FieldVariablePhiValue & phi() const override
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.
const MooseVariable * _temperature
Temperature, if provided. This is used only to get the trial functions.
std::vector< std::vector< const MaterialProperty< RankTwoTensor > * > > _deigenstrain_dargs
Eigenstrain derivatives wrt generate coupleds.
const MaterialProperty< RankTwoTensor > & _f_inv
The inverse increment deformation gradient.

◆ computeQpOffDiagJacobian()

Real LagrangianStressDivergenceBase::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtualinherited

Reimplemented in TotalLagrangianWeakPlaneStress.

Definition at line 121 of file LagrangianStressDivergenceBase.C.

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 unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
virtual Real computeQpJacobianDisplacement(unsigned int alpha, unsigned int beta)=0
unsigned int number() const
virtual Real computeQpJacobianOutOfPlaneStrain()=0
const MooseVariable * _out_of_plane_strain
Out-of-plane strain, if provided.
std::vector< unsigned int > _disp_nums
The displacement numbers.
const unsigned int _ndisp
Total number of displacements/size of residual vector.
const MooseVariable * _temperature
Temperature, if provided. This is used only to get the trial functions.
IntRange< T > make_range(T beg, T end)
virtual Real computeQpJacobianTemperature(unsigned int cvar)=0

◆ computeQpOffDiagJacobianScalar()

Real HomogenizedTotalLagrangianStressDivergence::computeQpOffDiagJacobianScalar ( const unsigned int  svar_num)
overrideprotectedvirtual

Method for computing d-_var-residual / d-svar at quadrature points.

Definition at line 217 of file HomogenizedTotalLagrangianStressDivergence.C.

Referenced by computeOffDiagJacobianScalarLocal().

219 {
220  return _dpk1[_qp].contractionKl(_m, _n, gradTest(_alpha));
221 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
virtual RankTwoTensor gradTest(unsigned int component) override
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.
unsigned int _m
Indices for off-diagonal Jacobian components.

◆ computeQpResidual()

template<class G >
Real TotalLagrangianStressDivergenceBase< G >::computeQpResidual ( )
overrideprotectedvirtualinherited

Reimplemented in TotalLagrangianWeakPlaneStress.

Definition at line 94 of file TotalLagrangianStressDivergenceBase.C.

95 {
96  return gradTest(_alpha).doubleContraction(_pk1[_qp]);
97 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
const MaterialProperty< RankTwoTensor > & _pk1
The 1st Piola-Kirchhoff stress.
virtual RankTwoTensor gradTest(unsigned int component) override
Real doubleContraction(const RankTwoTensorTempl< Real > &a) const

◆ computeScalarJacobian()

void HomogenizedTotalLagrangianStressDivergence::computeScalarJacobian ( )
overrideprotectedvirtual

Method for computing the scalar variable part of Jacobian for d-_kappa-residual / d-_kappa.

Definition at line 98 of file HomogenizedTotalLagrangianStressDivergence.C.

99 {
100  if (_alpha != 0)
101  return;
102 
103  _local_ke.resize(_k_order, _k_order);
104 
105  // only assemble scalar residual once; i.e. when handling the first displacement component
106  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
107  {
108  initScalarQpJacobian(_kappa_var);
109  const auto dV = _JxW[_qp] * _coord[_qp];
110 
111  // index for Jacobian row
112  unsigned int h = 0;
113 
114  for (const auto & [indices1, constraint1] : cmap())
115  {
116  const auto [i, j] = indices1;
117  const auto ctype = constraint1.first;
118 
119  // index for Jacobian col
120  unsigned int m = 0;
121 
122  for (const auto & [indices2, constraint2] : cmap())
123  {
124  const auto [k, l] = indices2;
126  _local_ke(h, m++) += dV * (_dpk1[_qp](i, j, k, l));
127  else if (ctype == Homogenization::ConstraintType::Strain)
128  {
129  if (_large_kinematics)
130  _local_ke(h, m++) += dV * (Real(i == k && j == l));
131  else
132  _local_ke(h, m++) += dV * (0.5 * Real(i == k && j == l) + 0.5 * Real(i == l && j == k));
133  }
134  else
135  mooseError("Unknown constraint type in Jacobian calculator!");
136  }
137  h++;
138  }
139  }
140 
141  addJacobian(_assembly,
142  _local_ke,
143  _kappa_var_ptr->dofIndices(),
144  _kappa_var_ptr->dofIndices(),
145  _kappa_var_ptr->scalingFactor());
146 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
void mooseError(Args &&... args)
const bool _large_kinematics
If true use large deformation kinematics.
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
static const std::string k
Definition: NS.h:130
const Homogenization::ConstraintMap & cmap() const
Get the constraint map.

◆ computeScalarOffDiagJacobian()

void HomogenizedTotalLagrangianStressDivergence::computeScalarOffDiagJacobian ( const unsigned int  jvar_num)
overrideprotectedvirtual

Method for computing an off-diagonal jacobian component d-_kappa-residual / d-jvar.

jvar is looped over all field variables, which herein is just disp_x and disp_y

Definition at line 149 of file HomogenizedTotalLagrangianStressDivergence.C.

151 {
152  // ONLY assemble the contribution from _alpha component, which is connected with _var
153  // The other components are handled by other kernel instances with other _alpha
154  if (jvar_num != _var.number())
155  return;
156 
157  const auto & jvar = getVariable(jvar_num);
158  const auto jvar_size = jvar.phiSize();
159  _local_ke.resize(_k_order, jvar_size);
160 
161  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
162  {
163  const auto dV = _JxW[_qp] * _coord[_qp];
164 
165  // index for Jacobian row
166  unsigned int h = 0;
167 
168  for (const auto & [indices, constraint] : cmap())
169  {
170  std::tie(_m, _n) = indices;
171  _ctype = constraint.first;
172  initScalarQpOffDiagJacobian(jvar);
173  for (_j = 0; _j < jvar_size; _j++)
174  _local_ke(h, _j) += dV * computeScalarQpOffDiagJacobian(jvar_num);
175  h++;
176  }
177  }
178 
179  addJacobian(_assembly,
180  _local_ke,
181  _kappa_var_ptr->dofIndices(),
182  jvar.dofIndices(),
183  _kappa_var_ptr->scalingFactor());
184 }
virtual Real computeScalarQpOffDiagJacobian(const unsigned int jvar_num) override
Method for computing an off-diagonal jacobian component at quadrature points.
unsigned int _m
Indices for off-diagonal Jacobian components.
Homogenization::ConstraintType _ctype
Type of current homogenization constraint.
const Homogenization::ConstraintMap & cmap() const
Get the constraint map.

◆ computeScalarQpOffDiagJacobian()

Real HomogenizedTotalLagrangianStressDivergence::computeScalarQpOffDiagJacobian ( const unsigned int  jvar_num)
overrideprotectedvirtual

Method for computing an off-diagonal jacobian component at quadrature points.

Definition at line 224 of file HomogenizedTotalLagrangianStressDivergence.C.

Referenced by computeScalarOffDiagJacobian().

226 {
228  return _dpk1[_qp].contractionIj(_m, _n, gradTrial(_alpha));
230  if (_large_kinematics)
231  return Real(_m == _alpha) * gradTrial(_alpha)(_m, _n);
232  else
233  return 0.5 * (Real(_m == _alpha) * gradTrial(_alpha)(_m, _n) +
234  Real(_n == _alpha) * gradTrial(_alpha)(_n, _m));
235  else
236  mooseError("Unknown constraint type in kernel calculation!");
237 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
void mooseError(Args &&... args)
const bool _large_kinematics
If true use large deformation kinematics.
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
unsigned int _m
Indices for off-diagonal Jacobian components.
Homogenization::ConstraintType _ctype
Type of current homogenization constraint.
virtual RankTwoTensor gradTrial(unsigned int component) override

◆ computeScalarResidual()

void HomogenizedTotalLagrangianStressDivergence::computeScalarResidual ( )
overrideprotectedvirtual

Method for computing the scalar part of residual for _kappa.

Definition at line 46 of file HomogenizedTotalLagrangianStressDivergence.C.

47 {
48  if (_alpha != 0)
49  return;
50 
51  std::vector<Real> scalar_residuals(_k_order);
52 
53  // only assemble scalar residual once; i.e. when handling the first displacement component
54  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
55  {
56  initScalarQpResidual();
57  const auto dV = _JxW[_qp] * _coord[_qp];
58 
59  // index for residual vector
60  unsigned int h = 0;
61 
62  for (const auto & [indices, constraint] : cmap())
63  {
64  const auto [i, j] = indices;
65  const auto [ctype, ctarget] = constraint;
66  const auto cval = ctarget->value(_t, _q_point[_qp]);
67 
68  // value to be constrained
69  Real val;
71  {
73  val = _pk1[_qp](i, j);
75  val = _F[_qp](i, j) - (Real(i == j));
76  else
77  mooseError("Unknown constraint type in the integral!");
78  }
79  else
80  {
82  val = _pk1[_qp](i, j);
84  val = 0.5 * (_F[_qp](i, j) + _F[_qp](j, i)) - (Real(i == j));
85  else
86  mooseError("Unknown constraint type in the integral!");
87  }
88 
89  scalar_residuals[h++] += (val - cval) * dV;
90  }
91  }
92 
93  addResiduals(
94  _assembly, scalar_residuals, _kappa_var_ptr->dofIndices(), _kappa_var_ptr->scalingFactor());
95 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
const MaterialProperty< RankTwoTensor > & _pk1
The 1st Piola-Kirchhoff stress.
void mooseError(Args &&... args)
const MaterialProperty< RankTwoTensor > & _F
The actual (stabilized) deformation gradient.
const bool _large_kinematics
If true use large deformation kinematics.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const Homogenization::ConstraintMap & cmap() const
Get the constraint map.

◆ gradTest()

template<class G >
RankTwoTensor TotalLagrangianStressDivergenceBase< G >::gradTest ( unsigned int  component)
overrideprotectedvirtualinherited

Implements LagrangianStressDivergenceBase.

Definition at line 33 of file TotalLagrangianStressDivergenceBase.C.

Referenced by computeQpOffDiagJacobianScalar().

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 }
static const std::string component
Definition: NS.h:153

◆ gradTrial()

template<class G >
RankTwoTensor TotalLagrangianStressDivergenceBase< G >::gradTrial ( unsigned int  component)
overrideprotectedvirtualinherited

Implements LagrangianStressDivergenceBase.

Definition at line 41 of file TotalLagrangianStressDivergenceBase.C.

Referenced by TotalLagrangianWeakPlaneStress::computeQpOffDiagJacobian(), and computeScalarQpOffDiagJacobian().

42 {
44 }
virtual RankTwoTensor gradTrialUnstabilized(unsigned int component)
The unstabilized trial function gradient.
static const std::string component
Definition: NS.h:153
const bool _stabilize_strain
If true calculate the deformation gradient derivatives for F_bar.
virtual RankTwoTensor gradTrialStabilized(unsigned int component)
The stabilized trial function gradient.

◆ initialSetup() [1/4]

template<>
void TotalLagrangianStressDivergenceBase< GradientOperatorCartesian >::initialSetup ( )
inlineinherited

Definition at line 26 of file TotalLagrangianStressDivergence.h.

27 {
28  if (getBlockCoordSystem() != Moose::COORD_XYZ)
29  mooseError("This kernel should only act in Cartesian coordinates.");
30 }
void mooseError(Args &&... args)

◆ initialSetup() [2/4]

Definition at line 26 of file TotalLagrangianStressDivergenceAxisymmetricCylindrical.h.

27 {
28  if (getBlockCoordSystem() != Moose::COORD_RZ)
29  mooseError("This kernel should only act in axisymmetric cylindrical coordinates.");
30 }
void mooseError(Args &&... args)

◆ initialSetup() [3/4]

Definition at line 26 of file TotalLagrangianStressDivergenceCentrosymmetricSpherical.h.

27 {
28  if (getBlockCoordSystem() != Moose::COORD_RSPHERICAL)
29  mooseError("This kernel should only act in centrosymmetric spherical coordinates.");
30 }
void mooseError(Args &&... args)
COORD_RSPHERICAL

◆ initialSetup() [4/4]

template<class G >
virtual void TotalLagrangianStressDivergenceBase< G >::initialSetup ( )
overridevirtualinherited

◆ precalculateJacobian()

void LagrangianStressDivergenceBase::precalculateJacobian ( )
overrideprotectedvirtualinherited

Reimplemented in TotalLagrangianWeakPlaneStress.

Definition at line 85 of file LagrangianStressDivergenceBase.C.

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 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
virtual void precalculateJacobianDisplacement(unsigned int component)=0
Prepare the average shape function gradients for stabilization.
const bool _stabilize_strain
If true calculate the deformation gradient derivatives for F_bar.

◆ precalculateJacobianDisplacement()

template<class G >
void TotalLagrangianStressDivergenceBase< G >::precalculateJacobianDisplacement ( unsigned int  component)
overrideprotectedvirtualinherited

Prepare the average shape function gradients for stabilization.

Implements LagrangianStressDivergenceBase.

Definition at line 80 of file TotalLagrangianStressDivergenceBase.C.

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 }
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
static const std::string component
Definition: NS.h:153
IntRange< T > make_range(T beg, T end)
auto elementAverage(const Functor &f, const MooseArray< Real > &JxW, const MooseArray< Real > &coord)
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")

◆ precalculateOffDiagJacobian()

void LagrangianStressDivergenceBase::precalculateOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtualinherited

Definition at line 98 of file LagrangianStressDivergenceBase.C.

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 }
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
virtual void precalculateJacobianDisplacement(unsigned int component)=0
Prepare the average shape function gradients for stabilization.
const bool _stabilize_strain
If true calculate the deformation gradient derivatives for F_bar.
std::vector< unsigned int > _disp_nums
The displacement numbers.
const unsigned int _ndisp
Total number of displacements/size of residual vector.
IntRange< T > make_range(T beg, T end)

◆ validParams()

InputParameters HomogenizedTotalLagrangianStressDivergence::validParams ( )
static

Definition at line 19 of file HomogenizedTotalLagrangianStressDivergence.C.

20 {
22  params.addClassDescription("Total Lagrangian stress equilibrium kernel with "
23  "homogenization constraint Jacobian terms");
24  params.renameCoupledVar(
25  "scalar_variable", "macro_var", "Optional scalar field with the macro gradient");
26 
27  return params;
28 }
void renameCoupledVar(const std::string &old_name, const std::string &new_name, const std::string &new_docstring)
static InputParameters validParams()
void addClassDescription(const std::string &doc_string)

Member Data Documentation

◆ _alpha

const unsigned int LagrangianStressDivergenceBase::_alpha
protectedinherited

◆ _avg_grad_trial

std::vector<std::vector<RankTwoTensor> > LagrangianStressDivergenceBase::_avg_grad_trial
protectedinherited

◆ _base_name

const std::string LagrangianStressDivergenceBase::_base_name
protectedinherited

Prepend to the material properties.

Definition at line 69 of file LagrangianStressDivergenceBase.h.

◆ _ctype

Homogenization::ConstraintType HomogenizedTotalLagrangianStressDivergence::_ctype = Homogenization::ConstraintType::None
private

Type of current homogenization constraint.

Definition at line 65 of file HomogenizedTotalLagrangianStressDivergence.h.

Referenced by computeOffDiagJacobianScalarLocal(), computeScalarOffDiagJacobian(), and computeScalarQpOffDiagJacobian().

◆ _deigenstrain_dargs

std::vector<std::vector<const MaterialProperty<RankTwoTensor> *> > LagrangianStressDivergenceBase::_deigenstrain_dargs
protectedinherited

Eigenstrain derivatives wrt generate coupleds.

Definition at line 107 of file LagrangianStressDivergenceBase.h.

Referenced by LagrangianStressDivergenceBase::LagrangianStressDivergenceBase().

◆ _disp_nums

std::vector<unsigned int> LagrangianStressDivergenceBase::_disp_nums
protectedinherited

◆ _dpk1

template<class G >
const MaterialProperty<RankFourTensor>& TotalLagrangianStressDivergenceBase< G >::_dpk1
protectedinherited

◆ _F

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_F
protectedinherited

The actual (stabilized) deformation gradient.

Definition at line 98 of file LagrangianStressDivergenceBase.h.

Referenced by computeScalarResidual().

◆ _F_avg

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_F_avg
protectedinherited

The element-average deformation gradient.

Definition at line 89 of file LagrangianStressDivergenceBase.h.

◆ _f_inv

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_f_inv
protectedinherited

The inverse increment deformation gradient.

Definition at line 92 of file LagrangianStressDivergenceBase.h.

◆ _F_inv

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_F_inv
protectedinherited

The inverse deformation gradient.

Definition at line 95 of file LagrangianStressDivergenceBase.h.

◆ _F_ust

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBase::_F_ust
protectedinherited

The unmodified deformation gradient.

Definition at line 86 of file LagrangianStressDivergenceBase.h.

◆ _large_kinematics

const bool LagrangianStressDivergenceBase::_large_kinematics
protectedinherited

◆ _m

unsigned int HomogenizedTotalLagrangianStressDivergence::_m
private

◆ _n

unsigned int HomogenizedTotalLagrangianStressDivergence::_n
private

◆ _ndisp

const unsigned int LagrangianStressDivergenceBase::_ndisp
protectedinherited

◆ _out_of_plane_strain

const MooseVariable* LagrangianStressDivergenceBase::_out_of_plane_strain
protectedinherited

◆ _pk1

template<class G >
const MaterialProperty<RankTwoTensor>& TotalLagrangianStressDivergenceBase< G >::_pk1
protectedinherited

The 1st Piola-Kirchhoff stress.

Definition at line 45 of file TotalLagrangianStressDivergenceBase.h.

Referenced by TotalLagrangianWeakPlaneStress::computeQpResidual(), and computeScalarResidual().

◆ _stabilize_strain

const bool LagrangianStressDivergenceBase::_stabilize_strain
protectedinherited

If true calculate the deformation gradient derivatives for F_bar.

Definition at line 66 of file LagrangianStressDivergenceBase.h.

Referenced by LagrangianStressDivergenceBase::precalculateJacobian(), and LagrangianStressDivergenceBase::precalculateOffDiagJacobian().

◆ _temperature

const MooseVariable* LagrangianStressDivergenceBase::_temperature
protectedinherited

Temperature, if provided. This is used only to get the trial functions.

Definition at line 101 of file LagrangianStressDivergenceBase.h.

Referenced by LagrangianStressDivergenceBase::computeQpOffDiagJacobian().


The documentation for this class was generated from the following files: