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

Total Lagrangian formulation with most homogenization terms (one disp_xyz field and one scalar) The macro_gradient variable is split into two scalars: the first component called '_hvar' herein and all other components called '_avar' herein. More...

#include <HomogenizedTotalLagrangianStressDivergenceR.h>

Inheritance diagram for HomogenizedTotalLagrangianStressDivergenceR:
[legend]

Public Types

typedef std::vector< intJvarMap
 

Public Member Functions

 HomogenizedTotalLagrangianStressDivergenceR (const InputParameters &parameters)
 
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 Real computeQpResidual () override
 
virtual Real computeQpJacobianDisplacement (unsigned int alpha, unsigned int beta) override
 
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 jvar is looped over all field variables, which herein is just disp_x and disp_y. 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) override
 Method for computing d-_var-residual / d-_svar at quadrature points. More...
 
virtual void computeScalarOffDiagJacobianScalar (const unsigned int svar_num) override
 Method for computing an off-diagonal jacobian component d-_kappa-residual / d-svar svar is looped over other scalar variables, which herein is just _kappa_other. 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 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 unsigned int _beta
 Which component of the scalar vector residual this constraint is responsible for. More...
 
const MooseVariableScalar *const _kappao_var_ptr
 (Pointer to) Scalar variable this kernel operates on More...
 
const unsigned int _kappao_var
 The unknown scalar variable ID. More...
 
const unsigned int _ko_order
 Order of the scalar variable, used in several places. More...
 
const VariableValue_kappa_other
 Reference to the current solution at the current quadrature point. More...
 
HomogenizationR::ConstraintMap _cmap
 Type of each constraint (stress or strain) for each component. More...
 
HomogenizationR::ConstraintType _ctype = HomogenizationR::ConstraintType::None
 The constraint type; initialize with 'none'. More...
 
unsigned int _m
 Used internally to iterate over each scalar component. More...
 
unsigned int _n
 
unsigned int _a
 
unsigned int _b
 
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
 

Detailed Description

Total Lagrangian formulation with most homogenization terms (one disp_xyz field and one scalar) The macro_gradient variable is split into two scalars: the first component called '_hvar' herein and all other components called '_avar' herein.

For parameter _beta = 0, the primary scalar (_kappa) is _hvar and the coupled scalar is _avar. For parameter _beta = 1, the primary scalar (_kappa) is _avar and the coupled scalar is _hvar. Just like the primary field variable (_var) is either disp_x or disp_y or disp_z depending on _alpha.

Thus, each instance of HomogenizedTotalLagrangianStressDivergenceR acts on one field variable (_disp_alpha) and one scalar variable (_hvar_beta). The job of the kernel is to assemble the residual of all dofs of _disp_alpha and of all dofs of _hvar_beta (namely, selected rows). Also, it assembles the ENTIRE row for _disp_alpha and _hvar_beta (namely the columns from all dofs of all _disp field variables and all dofs of all scalar variables _hvar and _avar). The rows for the other field/scalar variables are handled by other instances of the kernel, according to the flags compute_scalar_residuals and compute_field_residuals. When compute_field_residuals is given, only component=_alpha matters and beta = {0,1} is looped. When compute_scalar_residuals is given, only prime_scalar=_beta matters and alpha = {0,1,2} is looped.

In summary, for x=disp_x etc. and h=_hvar and a=_avar, then the contributions of the instances are _alpha=0 R = [Rx, 00, 00, 00, 00 ]^T J = [Jxx, Jxy, Jxz, Jxh, Jxa] _alpha=1 R = [00, Ry, 00, 00, 00 ]^T J = [Jyx, Jyy, Jyz, Jyh, Jya] _alpha=2 R = [00, 00, Rz, 00, 00 ]^T J = [Jzx, Jzy, Jzz, Jzh, Jza] _beta=0 R = [00, 00, 00, Rh, 00 ]^T J = [Jhx, Jhy, Jhz, Jhh, Jha] _beta=1 R = [00, 00, 00, 00, Ra ]^T J = [Jax, Jay, Jaz, Jah, Jaa]

In this manner, the full R and J are obtained with NO duplication of jobs: R = [Rx, Ry, Rz, Rh, Ra ]^T J = [Jxx, Jxy, Jxz, Jxh, Jxa Jyx, Jyy, Jyz, Jyh, Jya Jzx, Jzy, Jzz, Jzh, Jza Jhx, Jhy, Jhz, Jhh, Jha Jax, Jay, Jaz, Jah, Jaa]

Definition at line 72 of file HomogenizedTotalLagrangianStressDivergenceR.h.

Constructor & Destructor Documentation

◆ HomogenizedTotalLagrangianStressDivergenceR()

HomogenizedTotalLagrangianStressDivergenceR::HomogenizedTotalLagrangianStressDivergenceR ( const InputParameters parameters)

Definition at line 45 of file HomogenizedTotalLagrangianStressDivergenceR.C.

48  _beta(getParam<unsigned int>("prime_scalar")),
49  _kappao_var_ptr(getScalarVar("macro_other", 0)),
50  _kappao_var(coupledScalar("macro_other")),
51  _ko_order(getScalarVar("macro_other", 0)->order()),
52  _kappa_other(coupledScalarValue("macro_other"))
53 {
54  // Constraint types
55  auto types = getParam<MultiMooseEnum>("constraint_types");
56  if (types.size() != Moose::dim * Moose::dim)
57  mooseError("Number of constraint types must equal dim * dim. ", types.size(), " are provided.");
58 
59  // Targets to hit
60  const std::vector<FunctionName> & fnames = getParam<std::vector<FunctionName>>("targets");
61 
62  // Prepare the constraint map
63  unsigned int fcount = 0;
64  for (const auto j : make_range(Moose::dim))
65  for (const auto i : make_range(Moose::dim))
66  {
67  const auto idx = i + Moose::dim * j;
68  const auto ctype = static_cast<HomogenizationR::ConstraintType>(types.get(idx));
70  {
71  const Function * const f = &getFunctionByName(fnames[fcount++]);
72  _cmap[{i, j}] = {ctype, f};
73  }
74  }
75 }
void mooseError(Args &&... args)
static constexpr std::size_t dim
const unsigned int _kappao_var
The unknown scalar variable ID.
HomogenizationR::ConstraintMap _cmap
Type of each constraint (stress or strain) for each component.
Real f(Real x)
Test function for Brents method.
const unsigned int _ko_order
Order of the scalar variable, used in several places.
const MooseVariableScalar *const _kappao_var_ptr
(Pointer to) Scalar variable this kernel operates on
const VariableValue & _kappa_other
Reference to the current solution at the current quadrature point.
TotalLagrangianStressDivergenceBaseS< GradientOperatorCartesian > TotalLagrangianStressDivergenceS
IntRange< T > make_range(T beg, T end)
const unsigned int _beta
Which component of the scalar vector residual this constraint is responsible for. ...
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
ConstraintType
Constraint type: stress/PK stress or strain/deformation gradient.
unsigned int idx(const ElemType type, const unsigned int nx, const unsigned int i, const unsigned int j)

Member Function Documentation

◆ baseParams()

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

Definition at line 14 of file TotalLagrangianStressDivergenceBaseS.C.

Referenced by TotalLagrangianStressDivergenceBaseS< 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)

◆ computeOffDiagJacobianScalarLocal()

void HomogenizedTotalLagrangianStressDivergenceR::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 and _kappa_other

Definition at line 283 of file HomogenizedTotalLagrangianStressDivergenceR.C.

285 {
286  // Bail if jvar not coupled
287  if ((svar_num == _kappa_var) || (svar_num == _kappao_var))
288  {
289  // Get dofs and order of this scalar; at least one will be _kappa_var
290  const auto & svar = _sys.getScalarVariable(_tid, svar_num);
291  const unsigned int s_order = svar.order();
292  _local_ke.resize(_test.size(), s_order);
293 
294  // set the local beta based on the current svar_num
295  unsigned int beta = 0;
296  if (svar_num == _kappa_var)
297  beta = 0;
298  else // svar_num == _kappao_var
299  beta = 1;
300 
301  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
302  {
303  // single index for Jacobian row; double indices for constraint tensor component
304  unsigned int l = 0;
305  Real dV = _JxW[_qp] * _coord[_qp];
306  for (auto && [indices, constraint] : _cmap)
307  {
308  // identical logic to computeScalarResidual, but for column index
309  if (beta == 0)
310  {
311  if (l == 1)
312  break;
313  }
314  else
315  {
316  if (l == 0)
317  {
318  l++;
319  continue;
320  }
321  }
322  // copy constraint indices to protected variables to pass to Qp routine
323  _m = indices.first;
324  _n = indices.second;
325  _ctype = constraint.first;
326  initScalarQpJacobian(svar_num);
327  for (_i = 0; _i < _test.size(); _i++)
328  {
329  _local_ke(_i, -beta + l) += dV * computeQpOffDiagJacobianScalar(svar_num);
330  }
331  l++;
332  }
333  }
334 
335  addJacobian(_assembly, _local_ke, _var.dofIndices(), svar.dofIndices(), _var.scalingFactor());
336  }
337 }
const unsigned int _kappao_var
The unknown scalar variable ID.
HomogenizationR::ConstraintType _ctype
The constraint type; initialize with &#39;none&#39;.
HomogenizationR::ConstraintMap _cmap
Type of each constraint (stress or strain) for each component.
virtual Real computeQpOffDiagJacobianScalar(const unsigned int) override
Method for computing d-_var-residual / d-_svar at quadrature points.
unsigned int _m
Used internally to iterate over each scalar component.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real

◆ computeQpJacobian()

Real LagrangianStressDivergenceBaseS::computeQpJacobian ( )
overrideprotectedvirtualinherited

Definition at line 115 of file LagrangianStressDivergenceBaseS.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()

Real HomogenizedTotalLagrangianStressDivergenceR::computeQpJacobianDisplacement ( unsigned int  alpha,
unsigned int  beta 
)
overrideprotectedvirtual

Reimplemented from TotalLagrangianStressDivergenceBaseS< G >.

Definition at line 85 of file HomogenizedTotalLagrangianStressDivergenceR.C.

87 {
88  // Assemble J-alpha-beta
89  return gradTest(alpha).doubleContraction(_dpk1[_qp] * gradTrial(beta));
90 }
virtual RankTwoTensor gradTest(unsigned int component) override
Real doubleContraction(const RankTwoTensorTempl< Real > &a) const
static const std::string alpha
Definition: NS.h:134
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.
virtual RankTwoTensor gradTrial(unsigned int component) override

◆ computeQpJacobianOutOfPlaneStrain()

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

Implements LagrangianStressDivergenceBaseS.

Definition at line 128 of file TotalLagrangianStressDivergenceBaseS.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.
const FieldVariablePhiValue & phi() const override
virtual RankTwoTensor gradTest(unsigned int component) override
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.

◆ computeQpJacobianTemperature()

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

Implements LagrangianStressDivergenceBaseS.

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

◆ computeQpOffDiagJacobian()

Real LagrangianStressDivergenceBaseS::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtualinherited

Definition at line 121 of file LagrangianStressDivergenceBaseS.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 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.
unsigned int number() const
const MooseVariable * _out_of_plane_strain
Out-of-plane strain, if provided.
const unsigned int _ndisp
Total number of displacements/size of residual vector.
virtual Real computeQpJacobianOutOfPlaneStrain()=0
virtual Real computeQpJacobianDisplacement(unsigned int alpha, unsigned int beta)=0
virtual Real computeQpJacobianTemperature(unsigned int cvar)=0
IntRange< T > make_range(T beg, T end)

◆ computeQpOffDiagJacobianScalar()

Real HomogenizedTotalLagrangianStressDivergenceR::computeQpOffDiagJacobianScalar ( const unsigned  int)
overrideprotectedvirtual

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

Definition at line 340 of file HomogenizedTotalLagrangianStressDivergenceR.C.

Referenced by computeOffDiagJacobianScalarLocal().

342 {
343  return _dpk1[_qp].contractionKl(_m, _n, gradTest(_alpha));
344 }
const unsigned int _alpha
Which component of the vector residual this kernel is responsible for.
virtual RankTwoTensor gradTest(unsigned int component) override
unsigned int _m
Used internally to iterate over each scalar component.
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.

◆ computeQpResidual()

Real HomogenizedTotalLagrangianStressDivergenceR::computeQpResidual ( )
overrideprotectedvirtual

Reimplemented from TotalLagrangianStressDivergenceBaseS< G >.

Definition at line 78 of file HomogenizedTotalLagrangianStressDivergenceR.C.

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

◆ computeScalarJacobian()

void HomogenizedTotalLagrangianStressDivergenceR::computeScalarJacobian ( )
overrideprotectedvirtual

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

Definition at line 155 of file HomogenizedTotalLagrangianStressDivergenceR.C.

156 {
157  _local_ke.resize(_k_order, _k_order);
158 
159  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
160  {
161  initScalarQpJacobian(_kappa_var);
162  Real dV = _JxW[_qp] * _coord[_qp];
163 
164  _h = 0;
165  for (auto && [indices1, constraint1] : _cmap)
166  {
167  auto && [i, j] = indices1;
168  auto && ctype = constraint1.first;
169 
170  // identical logic to computeScalarResidual
171  if (_beta == 0)
172  {
173  if (_h == 1)
174  break;
175  }
176  else
177  {
178  if (_h == 0)
179  {
180  _h++;
181  continue;
182  }
183  }
184 
185  _l = 0;
186  for (auto && indices2 : _cmap)
187  {
188  auto && [a, b] = indices2.first;
189 
190  // identical logic to computeScalarResidual, but for column index
191  if (_beta == 0)
192  {
193  if (_l == 1)
194  break;
195  }
196  else
197  {
198  if (_l == 0)
199  {
200  _l++;
201  continue;
202  }
203  }
204 
205  unsigned int c_ind = -_beta + _l; // move 1 column left if _beta=1 for the other scalar
206  _l++; // increment the index before we forget
208  _local_ke(-_beta + _h, c_ind) += dV * (_dpk1[_qp](i, j, a, b));
209  else if (ctype == HomogenizationR::ConstraintType::Strain)
210  {
211  if (_large_kinematics)
212  _local_ke(-_beta + _h, c_ind) += dV * (Real(i == a && j == b));
213  else
214  _local_ke(-_beta + _h, c_ind) +=
215  dV * (0.5 * Real(i == a && j == b) + 0.5 * Real(i == b && j == a));
216  }
217  else
218  mooseError("Unknown constraint type in Jacobian calculator!");
219  }
220  _h++;
221  }
222  }
223 
224  addJacobian(_assembly,
225  _local_ke,
226  _kappa_var_ptr->dofIndices(),
227  _kappa_var_ptr->dofIndices(),
228  _kappa_var_ptr->scalingFactor());
229 }
const bool _large_kinematics
If true use large deformation kinematics.
void mooseError(Args &&... args)
HomogenizationR::ConstraintMap _cmap
Type of each constraint (stress or strain) for each component.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const unsigned int _beta
Which component of the scalar vector residual this constraint is responsible for. ...
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.

◆ computeScalarOffDiagJacobian()

void HomogenizedTotalLagrangianStressDivergenceR::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 232 of file HomogenizedTotalLagrangianStressDivergenceR.C.

234 {
235  // Bail if jvar not coupled
236  if (getJvarMap()[jvar_num] >= 0)
237  {
238 
239  const auto & jvar = getVariable(jvar_num);
240  _local_ke.resize(_k_order, _test.size());
241 
242  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
243  {
244  // single index for Jacobian column; double indices for constraint tensor component
245  unsigned int h = 0;
246  Real dV = _JxW[_qp] * _coord[_qp];
247  for (auto && [indices, constraint] : _cmap)
248  {
249  // identical logic to computeScalarResidual
250  if (_beta == 0)
251  {
252  if (h == 1)
253  break;
254  }
255  else
256  {
257  if (h == 0)
258  {
259  h++;
260  continue;
261  }
262  }
263  // copy constraint indices to protected variables to pass to Qp routine
264  _m = indices.first;
265  _n = indices.second;
266  _ctype = constraint.first;
267  initScalarQpOffDiagJacobian(_var);
268  for (_j = 0; _j < _test.size(); _j++)
269  _local_ke(-_beta + h, _j) += dV * computeScalarQpOffDiagJacobian(jvar_num);
270  h++;
271  }
272  }
273 
274  addJacobian(_assembly,
275  _local_ke,
276  _kappa_var_ptr->dofIndices(),
277  jvar.dofIndices(),
278  _kappa_var_ptr->scalingFactor());
279  }
280 }
HomogenizationR::ConstraintType _ctype
The constraint type; initialize with &#39;none&#39;.
HomogenizationR::ConstraintMap _cmap
Type of each constraint (stress or strain) for each component.
unsigned int _m
Used internally to iterate over each scalar component.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
virtual Real computeScalarQpOffDiagJacobian(const unsigned int jvar_num) override
Method for computing an off-diagonal jacobian component at quadrature points.
const unsigned int _beta
Which component of the scalar vector residual this constraint is responsible for. ...

◆ computeScalarOffDiagJacobianScalar()

void HomogenizedTotalLagrangianStressDivergenceR::computeScalarOffDiagJacobianScalar ( const unsigned int  svar_num)
overrideprotectedvirtual

Method for computing an off-diagonal jacobian component d-_kappa-residual / d-svar svar is looped over other scalar variables, which herein is just _kappa_other.

Definition at line 370 of file HomogenizedTotalLagrangianStressDivergenceR.C.

372 {
373  // Only do this for the other macro variable
374  if (svar_num == _kappao_var)
375  {
376  _local_ke.resize(_k_order, _ko_order);
377 
378  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
379  {
380  initScalarQpJacobian(_kappa_var);
381  Real dV = _JxW[_qp] * _coord[_qp];
382 
383  _h = 0;
384  for (auto && [indices1, constraint1] : _cmap)
385  {
386  auto && [i, j] = indices1;
387  auto && ctype = constraint1.first;
388 
389  // identical logic to computeScalarResidual
390  if (_beta == 0)
391  {
392  if (_h == 1)
393  break;
394  }
395  else
396  {
397  if (_h == 0)
398  {
399  _h++;
400  continue;
401  }
402  }
403 
404  _l = 0;
405  for (auto && indices2 : _cmap)
406  {
407  auto && [a, b] = indices2.first;
408 
409  // OPPOSITE logic/scalar from computeScalarResidual, AND for column index
410  if (_beta == 1)
411  {
412  if (_l == 1) // only assemble first=0 component of _hvar, then break the loop
413  break;
414  }
415  else
416  {
417  if (_l == 0) // skip first component (_hvar) & continue to "first" component of _avar
418  {
419  _l++;
420  continue;
421  }
422  }
423 
424  unsigned int c_ind =
425  -(1 - _beta) + _l; // DON'T move 1 column left if _beta=1 for the other scalar
426  _l++;
428  _local_ke(-_beta + _h, c_ind) += dV * (_dpk1[_qp](i, j, a, b));
429  else if (ctype == HomogenizationR::ConstraintType::Strain)
430  {
431  if (_large_kinematics)
432  _local_ke(-_beta + _h, c_ind) += dV * (Real(i == a && j == b));
433  else
434  _local_ke(-_beta + _h, c_ind) +=
435  dV * (0.5 * Real(i == a && j == b) + 0.5 * Real(i == b && j == a));
436  }
437  else
438  mooseError("Unknown constraint type in Jacobian calculator!");
439  }
440  _h++;
441  }
442  }
443 
444  addJacobian(_assembly,
445  _local_ke,
446  _kappa_var_ptr->dofIndices(),
448  _kappa_var_ptr->scalingFactor());
449  }
450 }
const bool _large_kinematics
If true use large deformation kinematics.
void mooseError(Args &&... args)
const unsigned int _kappao_var
The unknown scalar variable ID.
HomogenizationR::ConstraintMap _cmap
Type of each constraint (stress or strain) for each component.
const unsigned int _ko_order
Order of the scalar variable, used in several places.
const MooseVariableScalar *const _kappao_var_ptr
(Pointer to) Scalar variable this kernel operates on
virtual const std::vector< dof_id_type > & dofIndices() const
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const unsigned int _beta
Which component of the scalar vector residual this constraint is responsible for. ...
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.

◆ computeScalarQpOffDiagJacobian()

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

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

Definition at line 347 of file HomogenizedTotalLagrangianStressDivergenceR.C.

Referenced by computeScalarOffDiagJacobian().

348 {
349  // set the local alpha based on the current jvar_num
350  for (auto alpha : make_range(_ndisp))
351  {
352  if (jvar_num == _disp_nums[alpha])
353  {
355  return _dpk1[_qp].contractionIj(_m, _n, gradTrial(alpha));
357  if (_large_kinematics)
358  return Real(_m == alpha) * gradTrial(alpha)(_m, _n);
359  else
360  return 0.5 * (Real(_m == alpha) * gradTrial(alpha)(_m, _n) +
361  Real(_n == alpha) * gradTrial(alpha)(_n, _m));
362  else
363  mooseError("Unknown constraint type in kernel calculation!");
364  }
365  }
366  return 0.0;
367 }
std::vector< unsigned int > _disp_nums
The displacement numbers.
const bool _large_kinematics
If true use large deformation kinematics.
void mooseError(Args &&... args)
HomogenizationR::ConstraintType _ctype
The constraint type; initialize with &#39;none&#39;.
const unsigned int _ndisp
Total number of displacements/size of residual vector.
unsigned int _m
Used internally to iterate over each scalar component.
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)
const MaterialProperty< RankFourTensor > & _dpk1
The derivative of the PK1 stress with respect to the deformation gradient.
virtual RankTwoTensor gradTrial(unsigned int component) override

◆ computeScalarResidual()

void HomogenizedTotalLagrangianStressDivergenceR::computeScalarResidual ( )
overrideprotectedvirtual

Method for computing the scalar part of residual for _kappa.

Definition at line 93 of file HomogenizedTotalLagrangianStressDivergenceR.C.

94 {
95  std::vector<Real> scalar_residuals(_k_order);
96 
97  for (_qp = 0; _qp < _qrule->n_points(); _qp++)
98  {
99  initScalarQpResidual();
100  Real dV = _JxW[_qp] * _coord[_qp];
101  _h = 0; // single index for residual vector; double indices for constraint tensor component
102  for (auto && [indices, constraint] : _cmap)
103  {
104  auto && [i, j] = indices;
105  auto && [ctype, ctarget] = constraint;
106 
107  // ONLY the component(s) that this constraint will contribute here;
108  // other one is handled in the other constraint
109  if (_beta == 0)
110  {
111  if (_h == 1) // only assemble first=0 component of _hvar, then break the loop
112  break;
113  }
114  else
115  {
116  // skip the first component (_hvar) and continue to "first" component of _avar
117  if (_h == 0)
118  {
119  _h++;
120  continue;
121  }
122  }
123 
124  // I am not great with C++ precedence; so, store the index
125  unsigned int r_ind = -_beta + _h; // move 1 row up if _beta=1 for the other scalar
126  _h++; // increment the index before we forget
127  if (_large_kinematics)
128  {
130  scalar_residuals[r_ind] += dV * (_pk1[_qp](i, j) - ctarget->value(_t, _q_point[_qp]));
131  else if (ctype == HomogenizationR::ConstraintType::Strain)
132  scalar_residuals[r_ind] +=
133  dV * (_F[_qp](i, j) - (Real(i == j) + ctarget->value(_t, _q_point[_qp])));
134  else
135  mooseError("Unknown constraint type in the integral!");
136  }
137  else
138  {
140  scalar_residuals[r_ind] += dV * (_pk1[_qp](i, j) - ctarget->value(_t, _q_point[_qp]));
141  else if (ctype == HomogenizationR::ConstraintType::Strain)
142  scalar_residuals[r_ind] += dV * (0.5 * (_F[_qp](i, j) + _F[_qp](j, i)) -
143  (Real(i == j) + ctarget->value(_t, _q_point[_qp])));
144  else
145  mooseError("Unknown constraint type in the integral!");
146  }
147  }
148  }
149 
150  addResiduals(
151  _assembly, scalar_residuals, _kappa_var_ptr->dofIndices(), _kappa_var_ptr->scalingFactor());
152 }
const bool _large_kinematics
If true use large deformation kinematics.
void mooseError(Args &&... args)
HomogenizationR::ConstraintMap _cmap
Type of each constraint (stress or strain) for each component.
const MaterialProperty< RankTwoTensor > & _F
The actual (stabilized) deformation gradient.
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
const unsigned int _beta
Which component of the scalar vector residual this constraint is responsible for. ...
static const std::complex< double > j(0, 1)
Complex number "j" (also known as "i")
const MaterialProperty< RankTwoTensor > & _pk1
The 1st Piola-Kirchhoff stress.

◆ gradTest()

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

◆ gradTrial()

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

Implements LagrangianStressDivergenceBaseS.

Definition at line 41 of file TotalLagrangianStressDivergenceBaseS.C.

Referenced by computeQpJacobianDisplacement(), HomogenizedTotalLagrangianStressDivergenceA::computeQpJacobianDisplacement(), HomogenizedTotalLagrangianStressDivergenceS::computeScalarQpOffDiagJacobian(), computeScalarQpOffDiagJacobian(), and HomogenizedTotalLagrangianStressDivergenceA::computeScalarQpOffDiagJacobian().

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

◆ initialSetup() [1/4]

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

Definition at line 26 of file TotalLagrangianStressDivergenceS.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 TotalLagrangianStressDivergenceCentrosymmetricSphericalS.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() [3/4]

Definition at line 26 of file TotalLagrangianStressDivergenceAxisymmetricCylindricalS.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() [4/4]

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

◆ precalculateJacobian()

void LagrangianStressDivergenceBaseS::precalculateJacobian ( )
overrideprotectedvirtualinherited

Definition at line 85 of file LagrangianStressDivergenceBaseS.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.
const bool _stabilize_strain
If true calculate the deformation gradient derivatives for F_bar.
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
virtual void precalculateJacobianDisplacement(unsigned int component)=0
Prepare the average shape function gradients for stabilization.

◆ precalculateJacobianDisplacement()

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

Prepare the average shape function gradients for stabilization.

Implements LagrangianStressDivergenceBaseS.

Definition at line 80 of file TotalLagrangianStressDivergenceBaseS.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 }
static const std::string component
Definition: NS.h:153
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
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 LagrangianStressDivergenceBaseS::precalculateOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtualinherited

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

◆ validParams()

InputParameters HomogenizedTotalLagrangianStressDivergenceR::validParams ( )
static

Definition at line 25 of file HomogenizedTotalLagrangianStressDivergenceR.C.

26 {
28  params.addClassDescription("Total Lagrangian stress equilibrium kernel with "
29  "homogenization constraint Jacobian terms");
30  params.renameCoupledVar(
31  "scalar_variable", "macro_var", "Optional scalar field with the macro gradient");
32  params.addRequiredCoupledVar("macro_other", "Other components of coupled scalar variable");
33  params.addRequiredParam<unsigned int>("prime_scalar", "Either 0=_var or 1=_other scalar");
35  "constraint_types",
37  "Type of each constraint: strain, stress, or none. The types are specified in the "
38  "column-major order, and there must be 9 entries in total.");
39  params.addRequiredParam<std::vector<FunctionName>>(
40  "targets", "Functions giving the targets to hit for constraint types that are not none.");
41 
42  return params;
43 }
static InputParameters validParams()
void renameCoupledVar(const std::string &old_name, const std::string &new_name, const std::string &new_docstring)
void addRequiredParam(const std::string &name, const std::string &doc_string)
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
void addClassDescription(const std::string &doc_string)
const MultiMooseEnum constraintType("strain stress none")
Moose constraint type, for input.

Member Data Documentation

◆ _a

unsigned int HomogenizedTotalLagrangianStressDivergenceR::_a
protected

◆ _alpha

const unsigned int LagrangianStressDivergenceBaseS::_alpha
protectedinherited

◆ _avg_grad_trial

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

◆ _b

unsigned int HomogenizedTotalLagrangianStressDivergenceR::_b
protected

◆ _base_name

const std::string LagrangianStressDivergenceBaseS::_base_name
protectedinherited

Prepend to the material properties.

Definition at line 69 of file LagrangianStressDivergenceBaseS.h.

◆ _beta

const unsigned int HomogenizedTotalLagrangianStressDivergenceR::_beta
protected

Which component of the scalar vector residual this constraint is responsible for.

Definition at line 123 of file HomogenizedTotalLagrangianStressDivergenceR.h.

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

◆ _cmap

HomogenizationR::ConstraintMap HomogenizedTotalLagrangianStressDivergenceR::_cmap
protected

◆ _ctype

HomogenizationR::ConstraintType HomogenizedTotalLagrangianStressDivergenceR::_ctype = HomogenizationR::ConstraintType::None
protected

The constraint type; initialize with 'none'.

Definition at line 141 of file HomogenizedTotalLagrangianStressDivergenceR.h.

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

◆ _deigenstrain_dargs

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

Eigenstrain derivatives wrt generate coupleds.

Definition at line 107 of file LagrangianStressDivergenceBaseS.h.

Referenced by LagrangianStressDivergenceBaseS::LagrangianStressDivergenceBaseS().

◆ _disp_nums

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

◆ _dpk1

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

◆ _F

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_F
protectedinherited

◆ _F_avg

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_F_avg
protectedinherited

The element-average deformation gradient.

Definition at line 89 of file LagrangianStressDivergenceBaseS.h.

◆ _f_inv

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_f_inv
protectedinherited

The inverse increment deformation gradient.

Definition at line 92 of file LagrangianStressDivergenceBaseS.h.

◆ _F_inv

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_F_inv
protectedinherited

The inverse deformation gradient.

Definition at line 95 of file LagrangianStressDivergenceBaseS.h.

◆ _F_ust

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_F_ust
protectedinherited

The unmodified deformation gradient.

Definition at line 86 of file LagrangianStressDivergenceBaseS.h.

◆ _kappa_other

const VariableValue& HomogenizedTotalLagrangianStressDivergenceR::_kappa_other
protected

Reference to the current solution at the current quadrature point.

Definition at line 135 of file HomogenizedTotalLagrangianStressDivergenceR.h.

◆ _kappao_var

const unsigned int HomogenizedTotalLagrangianStressDivergenceR::_kappao_var
protected

The unknown scalar variable ID.

Definition at line 129 of file HomogenizedTotalLagrangianStressDivergenceR.h.

Referenced by computeOffDiagJacobianScalarLocal(), and computeScalarOffDiagJacobianScalar().

◆ _kappao_var_ptr

const MooseVariableScalar* const HomogenizedTotalLagrangianStressDivergenceR::_kappao_var_ptr
protected

(Pointer to) Scalar variable this kernel operates on

Definition at line 126 of file HomogenizedTotalLagrangianStressDivergenceR.h.

Referenced by computeScalarOffDiagJacobianScalar().

◆ _ko_order

const unsigned int HomogenizedTotalLagrangianStressDivergenceR::_ko_order
protected

Order of the scalar variable, used in several places.

Definition at line 132 of file HomogenizedTotalLagrangianStressDivergenceR.h.

Referenced by computeScalarOffDiagJacobianScalar().

◆ _large_kinematics

const bool LagrangianStressDivergenceBaseS::_large_kinematics
protectedinherited

◆ _m

unsigned int HomogenizedTotalLagrangianStressDivergenceR::_m
protected

◆ _n

unsigned int HomogenizedTotalLagrangianStressDivergenceR::_n
protected

◆ _ndisp

const unsigned int LagrangianStressDivergenceBaseS::_ndisp
protectedinherited

◆ _out_of_plane_strain

const MooseVariable* LagrangianStressDivergenceBaseS::_out_of_plane_strain
protectedinherited

Out-of-plane strain, if provided.

Definition at line 104 of file LagrangianStressDivergenceBaseS.h.

Referenced by LagrangianStressDivergenceBaseS::computeQpOffDiagJacobian().

◆ _pk1

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

◆ _stabilize_strain

const bool LagrangianStressDivergenceBaseS::_stabilize_strain
protectedinherited

If true calculate the deformation gradient derivatives for F_bar.

Definition at line 66 of file LagrangianStressDivergenceBaseS.h.

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

◆ _temperature

const MooseVariable* LagrangianStressDivergenceBaseS::_temperature
protectedinherited

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

Definition at line 101 of file LagrangianStressDivergenceBaseS.h.

Referenced by LagrangianStressDivergenceBaseS::computeQpOffDiagJacobian().


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