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

Base class of the "Lagrangian" kernel system. More...

#include <LagrangianStressDivergenceBaseS.h>

Inheritance diagram for LagrangianStressDivergenceBaseS:
[legend]

Public Types

typedef std::vector< intJvarMap
 

Public Member Functions

 LagrangianStressDivergenceBaseS (const InputParameters &parameters)
 
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 ()
 

Protected Member Functions

virtual RankTwoTensor gradTest (unsigned int component)=0
 
virtual RankTwoTensor gradTrial (unsigned int component)=0
 
virtual void precalculateJacobian () override
 
virtual void precalculateOffDiagJacobian (unsigned int jvar) override
 
virtual void precalculateJacobianDisplacement (unsigned int component)=0
 Prepare the average shape function gradients for stabilization. More...
 
virtual Real computeQpJacobian () override
 
virtual Real computeQpOffDiagJacobian (unsigned int jvar) override
 
virtual Real computeQpJacobianDisplacement (unsigned int alpha, unsigned int beta)=0
 
virtual Real computeQpJacobianTemperature (unsigned int cvar)=0
 
virtual Real computeQpJacobianOutOfPlaneStrain ()=0
 

Protected Attributes

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

Base class of the "Lagrangian" kernel system.

This class provides a common structure for the "new" tensor_mechanics kernel system. The goals for this new system are 1) Always-correct jacobians 2) A cleaner material interface

This class provides common input properties and helper methods, most of the math has to be done in the subclasses

Definition at line 27 of file LagrangianStressDivergenceBaseS.h.

Constructor & Destructor Documentation

◆ LagrangianStressDivergenceBaseS()

LagrangianStressDivergenceBaseS::LagrangianStressDivergenceBaseS ( const InputParameters parameters)

Definition at line 42 of file LagrangianStressDivergenceBaseS.C.

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")),
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 }
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.
const bool _large_kinematics
If true use large deformation kinematics.
void mooseError(Args &&... args)
const MaterialProperty< RankTwoTensor > & _f_inv
The inverse increment deformation gradient.
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
const MaterialProperty< RankTwoTensor > & _F
The actual (stabilized) deformation gradient.
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.
const MaterialProperty< RankTwoTensor > & _F_inv
The inverse deformation gradient.
std::vector< std::vector< RankTwoTensor > > _avg_grad_trial
const std::string _base_name
Prepend to the material properties.
const MaterialProperty< RankTwoTensor > & _F_ust
The unmodified deformation gradient.
const MaterialProperty< RankTwoTensor > & _F_avg
The element-average deformation gradient.

Member Function Documentation

◆ computeQpJacobian()

Real LagrangianStressDivergenceBaseS::computeQpJacobian ( )
overrideprotectedvirtual

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()

virtual Real LagrangianStressDivergenceBaseS::computeQpJacobianDisplacement ( unsigned int  alpha,
unsigned int  beta 
)
protectedpure virtual

◆ computeQpJacobianOutOfPlaneStrain()

virtual Real LagrangianStressDivergenceBaseS::computeQpJacobianOutOfPlaneStrain ( )
protectedpure virtual

◆ computeQpJacobianTemperature()

virtual Real LagrangianStressDivergenceBaseS::computeQpJacobianTemperature ( unsigned int  cvar)
protectedpure virtual

◆ computeQpOffDiagJacobian()

Real LagrangianStressDivergenceBaseS::computeQpOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtual

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)

◆ gradTest()

virtual RankTwoTensor LagrangianStressDivergenceBaseS::gradTest ( unsigned int  component)
protectedpure virtual

◆ gradTrial()

virtual RankTwoTensor LagrangianStressDivergenceBaseS::gradTrial ( unsigned int  component)
protectedpure virtual

◆ precalculateJacobian()

void LagrangianStressDivergenceBaseS::precalculateJacobian ( )
overrideprotectedvirtual

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()

virtual void LagrangianStressDivergenceBaseS::precalculateJacobianDisplacement ( unsigned int  component)
protectedpure virtual

Prepare the average shape function gradients for stabilization.

Implemented in TotalLagrangianStressDivergenceBaseS< G >.

Referenced by precalculateJacobian(), and precalculateOffDiagJacobian().

◆ precalculateOffDiagJacobian()

void LagrangianStressDivergenceBaseS::precalculateOffDiagJacobian ( unsigned int  jvar)
overrideprotectedvirtual

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 LagrangianStressDivergenceBaseS::validParams ( )
static

Definition at line 13 of file LagrangianStressDivergenceBaseS.C.

Referenced by TotalLagrangianStressDivergenceBaseS< G >::baseParams().

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 }
void addParam(const std::string &name, const std::initializer_list< typename T::value_type > &value, const std::string &doc_string)
void addRequiredParam(const std::string &name, const std::string &doc_string)
void addCoupledVar(const std::string &name, const std::string &doc_string)
void addRequiredCoupledVar(const std::string &name, const std::string &doc_string)
static InputParameters validParams()

Member Data Documentation

◆ _alpha

const unsigned int LagrangianStressDivergenceBaseS::_alpha
protected

◆ _avg_grad_trial

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

◆ _base_name

const std::string LagrangianStressDivergenceBaseS::_base_name
protected

Prepend to the material properties.

Definition at line 69 of file LagrangianStressDivergenceBaseS.h.

◆ _deigenstrain_dargs

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

Eigenstrain derivatives wrt generate coupleds.

Definition at line 107 of file LagrangianStressDivergenceBaseS.h.

Referenced by LagrangianStressDivergenceBaseS().

◆ _disp_nums

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

◆ _F

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_F
protected

◆ _F_avg

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_F_avg
protected

The element-average deformation gradient.

Definition at line 89 of file LagrangianStressDivergenceBaseS.h.

◆ _f_inv

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_f_inv
protected

The inverse increment deformation gradient.

Definition at line 92 of file LagrangianStressDivergenceBaseS.h.

◆ _F_inv

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_F_inv
protected

The inverse deformation gradient.

Definition at line 95 of file LagrangianStressDivergenceBaseS.h.

◆ _F_ust

const MaterialProperty<RankTwoTensor>& LagrangianStressDivergenceBaseS::_F_ust
protected

The unmodified deformation gradient.

Definition at line 86 of file LagrangianStressDivergenceBaseS.h.

◆ _large_kinematics

const bool LagrangianStressDivergenceBaseS::_large_kinematics
protected

◆ _ndisp

const unsigned int LagrangianStressDivergenceBaseS::_ndisp
protected

◆ _out_of_plane_strain

const MooseVariable* LagrangianStressDivergenceBaseS::_out_of_plane_strain
protected

Out-of-plane strain, if provided.

Definition at line 104 of file LagrangianStressDivergenceBaseS.h.

Referenced by computeQpOffDiagJacobian().

◆ _stabilize_strain

const bool LagrangianStressDivergenceBaseS::_stabilize_strain
protected

If true calculate the deformation gradient derivatives for F_bar.

Definition at line 66 of file LagrangianStressDivergenceBaseS.h.

Referenced by precalculateJacobian(), and precalculateOffDiagJacobian().

◆ _temperature

const MooseVariable* LagrangianStressDivergenceBaseS::_temperature
protected

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

Definition at line 101 of file LagrangianStressDivergenceBaseS.h.

Referenced by computeQpOffDiagJacobian().


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