www.mooseframework.org
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Attributes | List of all members
ComputeFiniteStrain Class Reference

ComputeFiniteStrain defines a strain increment and rotation increment, for finite strains. More...

#include <ComputeFiniteStrain.h>

Inheritance diagram for ComputeFiniteStrain:
[legend]

Public Member Functions

 ComputeFiniteStrain (const InputParameters &parameters)
 
void computeProperties () override
 
void initialSetup () override
 

Static Public Member Functions

static InputParameters validParams ()
 
static MooseEnum decompositionType ()
 

Protected Member Functions

virtual void computeQpStrain ()
 
virtual void computeQpIncrements (RankTwoTensor &e, RankTwoTensor &r)
 
virtual void initQpStatefulProperties () override
 
void subtractEigenstrainIncrementFromStrain (RankTwoTensor &strain)
 
virtual void displacementIntegrityCheck ()
 

Protected Attributes

std::vector< RankTwoTensor_Fhat
 
std::vector< const VariableGradient * > _grad_disp_old
 
MaterialProperty< RankTwoTensor > & _strain_rate
 
MaterialProperty< RankTwoTensor > & _strain_increment
 
MaterialProperty< RankTwoTensor > & _rotation_increment
 
MaterialProperty< RankTwoTensor > & _deformation_gradient
 
const MaterialProperty< RankTwoTensor > & _mechanical_strain_old
 
const MaterialProperty< RankTwoTensor > & _total_strain_old
 
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains_old
 
unsigned int _ndisp
 Coupled displacement variables. More...
 
std::vector< const VariableValue * > _disp
 
std::vector< const VariableGradient * > _grad_disp
 
const std::string _base_name
 
MaterialProperty< RankTwoTensor > & _mechanical_strain
 
MaterialProperty< RankTwoTensor > & _total_strain
 
std::vector< MaterialPropertyName > _eigenstrain_names
 
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
 
const MaterialProperty< RankTwoTensor > * _global_strain
 
const bool _volumetric_locking_correction
 
const Real & _current_elem_volume
 

Private Types

enum  DecompMethod { DecompMethod::TaylorExpansion, DecompMethod::EigenSolution }
 

Private Attributes

const DecompMethod _decomposition_method
 

Detailed Description

ComputeFiniteStrain defines a strain increment and rotation increment, for finite strains.

Definition at line 22 of file ComputeFiniteStrain.h.

Member Enumeration Documentation

◆ DecompMethod

Enumerator
TaylorExpansion 
EigenSolution 

Definition at line 40 of file ComputeFiniteStrain.h.

41  {
42  TaylorExpansion,
43  EigenSolution
44  };

Constructor & Destructor Documentation

◆ ComputeFiniteStrain()

ComputeFiniteStrain::ComputeFiniteStrain ( const InputParameters &  parameters)

Definition at line 38 of file ComputeFiniteStrain.C.

39  : ComputeIncrementalStrainBase(parameters),
40  _Fhat(_fe_problem.getMaxQps()),
41  _decomposition_method(getParam<MooseEnum>("decomposition_method").getEnum<DecompMethod>())
42 {
43 }

Member Function Documentation

◆ computeProperties()

void ComputeFiniteStrain::computeProperties ( )
override

Definition at line 46 of file ComputeFiniteStrain.C.

47 {
48  RankTwoTensor ave_Fhat;
49  Real ave_dfgrd_det = 0.0;
50  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
51  {
52  // Deformation gradient
53  RankTwoTensor A((*_grad_disp[0])[_qp],
54  (*_grad_disp[1])[_qp],
55  (*_grad_disp[2])[_qp]); // Deformation gradient
56  RankTwoTensor Fbar((*_grad_disp_old[0])[_qp],
57  (*_grad_disp_old[1])[_qp],
58  (*_grad_disp_old[2])[_qp]); // Old Deformation gradient
59 
60  _deformation_gradient[_qp] = A;
61  _deformation_gradient[_qp].addIa(1.0); // Gauss point deformation gradient
62 
63  // A = gradU - gradUold
64  A -= Fbar;
65 
66  // Fbar = ( I + gradUold)
67  Fbar.addIa(1.0);
68 
69  // Incremental deformation gradient _Fhat = I + A Fbar^-1
70  _Fhat[_qp] = A * Fbar.inverse();
71  _Fhat[_qp].addIa(1.0);
72 
74  {
75  // Calculate average _Fhat for volumetric locking correction
76  ave_Fhat += _Fhat[_qp] * _JxW[_qp] * _coord[_qp];
77 
78  // Average deformation gradient
79  ave_dfgrd_det += _deformation_gradient[_qp].det() * _JxW[_qp] * _coord[_qp];
80  }
81  }
82 
84  {
85  // needed for volumetric locking correction
86  ave_Fhat /= _current_elem_volume;
87  // average deformation gradient
88  ave_dfgrd_det /= _current_elem_volume;
89  }
90 
91  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
92  {
94  {
95  // Finalize volumetric locking correction
96  _Fhat[_qp] *= std::cbrt(ave_Fhat.det() / _Fhat[_qp].det());
97  _deformation_gradient[_qp] *= std::cbrt(ave_dfgrd_det / _deformation_gradient[_qp].det());
98  }
99 
100  computeQpStrain();
101  }
102 }

◆ computeQpIncrements()

void ComputeFiniteStrain::computeQpIncrements ( RankTwoTensor e,
RankTwoTensor r 
)
protectedvirtual

Definition at line 137 of file ComputeFiniteStrain.C.

139 {
140  switch (_decomposition_method)
141  {
143  {
144  // inverse of _Fhat
145  const RankTwoTensor invFhat = _Fhat[_qp].inverse();
146 
147  // A = I - _Fhat^-1
148  RankTwoTensor A(RankTwoTensor::initIdentity);
149  A -= invFhat;
150 
151  // Cinv - I = A A^T - A - A^T;
152  RankTwoTensor Cinv_I = A * A.transpose() - A - A.transpose();
153 
154  // strain rate D from Taylor expansion, Chat = (-1/2(Chat^-1 - I) + 1/4*(Chat^-1 - I)^2 + ...
155  total_strain_increment = -Cinv_I * 0.5 + Cinv_I * Cinv_I * 0.25;
156 
157  const Real a[3] = {invFhat(1, 2) - invFhat(2, 1),
158  invFhat(2, 0) - invFhat(0, 2),
159  invFhat(0, 1) - invFhat(1, 0)};
160 
161  Real q = (a[0] * a[0] + a[1] * a[1] + a[2] * a[2]) / 4.0;
162  Real trFhatinv_1 = invFhat.trace() - 1.0;
163  const Real p = trFhatinv_1 * trFhatinv_1 / 4.0;
164 
165  // cos theta_a
166  const Real C1_squared = p +
167  3.0 * Utility::pow<2>(p) * (1.0 - (p + q)) / Utility::pow<2>(p + q) -
168  2.0 * Utility::pow<3>(p) * (1.0 - (p + q)) / Utility::pow<3>(p + q);
169  mooseAssert(C1_squared >= 0.0,
170  "Cannot take square root of a negative number. This may happen when elements "
171  "become heavily distorted.");
172  const Real C1 = std::sqrt(C1_squared);
173 
174  Real C2;
175  if (q > 0.01)
176  // (1-cos theta_a)/4q
177  C2 = (1.0 - C1) / (4.0 * q);
178  else
179  // alternate form for small q
180  C2 = 0.125 + q * 0.03125 * (Utility::pow<2>(p) - 12.0 * (p - 1.0)) / Utility::pow<2>(p) +
181  Utility::pow<2>(q) * (p - 2.0) * (Utility::pow<2>(p) - 10.0 * p + 32.0) /
182  Utility::pow<3>(p) +
183  Utility::pow<3>(q) *
184  (1104.0 - 992.0 * p + 376.0 * Utility::pow<2>(p) - 72.0 * Utility::pow<3>(p) +
185  5.0 * Utility::pow<4>(p)) /
186  (512.0 * Utility::pow<4>(p));
187 
188  const Real C3_test =
189  (p * q * (3.0 - q) + Utility::pow<3>(p) + Utility::pow<2>(q)) / Utility::pow<3>(p + q);
190  mooseAssert(C3_test >= 0.0,
191  "Cannot take square root of a negative number. This may happen when elements "
192  "become heavily distorted.");
193  const Real C3 = 0.5 * std::sqrt(C3_test); // sin theta_a/(2 sqrt(q))
194 
195  // Calculate incremental rotation. Note that this value is the transpose of that from Rashid,
196  // 93, so we transpose it before storing
197  RankTwoTensor R_incr;
198  R_incr.addIa(C1);
199  for (unsigned int i = 0; i < 3; ++i)
200  for (unsigned int j = 0; j < 3; ++j)
201  R_incr(i, j) += C2 * a[i] * a[j];
202 
203  R_incr(0, 1) += C3 * a[2];
204  R_incr(0, 2) -= C3 * a[1];
205  R_incr(1, 0) -= C3 * a[2];
206  R_incr(1, 2) += C3 * a[0];
207  R_incr(2, 0) += C3 * a[1];
208  R_incr(2, 1) -= C3 * a[0];
209 
210  rotation_increment = R_incr.transpose();
211  break;
212  }
213 
215  {
216  std::vector<Real> e_value(3);
217  RankTwoTensor e_vector, N1, N2, N3;
218 
219  RankTwoTensor Chat = _Fhat[_qp].transpose() * _Fhat[_qp];
220  Chat.symmetricEigenvaluesEigenvectors(e_value, e_vector);
221 
222  const Real lambda1 = std::sqrt(e_value[0]);
223  const Real lambda2 = std::sqrt(e_value[1]);
224  const Real lambda3 = std::sqrt(e_value[2]);
225 
226  N1.vectorOuterProduct(e_vector.column(0), e_vector.column(0));
227  N2.vectorOuterProduct(e_vector.column(1), e_vector.column(1));
228  N3.vectorOuterProduct(e_vector.column(2), e_vector.column(2));
229 
230  const RankTwoTensor Uhat = N1 * lambda1 + N2 * lambda2 + N3 * lambda3;
231  const RankTwoTensor invUhat(Uhat.inverse());
232 
233  rotation_increment = _Fhat[_qp] * invUhat;
234 
235  total_strain_increment =
236  N1 * std::log(lambda1) + N2 * std::log(lambda2) + N3 * std::log(lambda3);
237  break;
238  }
239 
240  default:
241  mooseError("ComputeFiniteStrain Error: Pass valid decomposition type: TaylorExpansion or "
242  "EigenSolution.");
243  }
244 }

Referenced by computeQpStrain().

◆ computeQpStrain()

void ComputeFiniteStrain::computeQpStrain ( )
protectedvirtual

Definition at line 105 of file ComputeFiniteStrain.C.

106 {
107  RankTwoTensor total_strain_increment;
108 
109  // two ways to calculate these increments: TaylorExpansion(default) or EigenSolution
110  computeQpIncrements(total_strain_increment, _rotation_increment[_qp]);
111 
112  _strain_increment[_qp] = total_strain_increment;
113 
114  // Remove the eigenstrain increment
116 
117  if (_dt > 0)
118  _strain_rate[_qp] = _strain_increment[_qp] / _dt;
119  else
120  _strain_rate[_qp].zero();
121 
122  // Update strain in intermediate configuration
124  _total_strain[_qp] = _total_strain_old[_qp] + total_strain_increment;
125 
126  // Rotate strain to current configuration
127  _mechanical_strain[_qp] =
128  _rotation_increment[_qp] * _mechanical_strain[_qp] * _rotation_increment[_qp].transpose();
129  _total_strain[_qp] =
130  _rotation_increment[_qp] * _total_strain[_qp] * _rotation_increment[_qp].transpose();
131 
132  if (_global_strain)
133  _total_strain[_qp] += (*_global_strain)[_qp];
134 }

Referenced by computeProperties(), Compute1DFiniteStrain::computeProperties(), Compute2DFiniteStrain::computeProperties(), and ComputeRSphericalFiniteStrain::computeProperties().

◆ decompositionType()

MooseEnum ComputeFiniteStrain::decompositionType ( )
static

Definition at line 17 of file ComputeFiniteStrain.C.

18 {
19  return MooseEnum("TaylorExpansion EigenSolution", "TaylorExpansion");
20 }

Referenced by TensorMechanicsActionBase::validParams(), and validParams().

◆ displacementIntegrityCheck()

void ComputeStrainBase::displacementIntegrityCheck ( )
protectedvirtualinherited

Reimplemented in Compute2DFiniteStrain, Compute2DIncrementalStrain, and Compute2DSmallStrain.

Definition at line 94 of file ComputeStrainBase.C.

95 {
96  // Checking for consistency between mesh size and length of the provided displacements vector
97  if (_ndisp != _mesh.dimension())
98  paramError(
99  "displacements",
100  "The number of variables supplied in 'displacements' must match the mesh dimension.");
101 }

Referenced by ComputeStrainBase::initialSetup().

◆ initialSetup()

void ComputeIncrementalStrainBase::initialSetup ( )
overrideinherited

Definition at line 38 of file ComputeIncrementalStrainBase.C.

39 {
41  for (unsigned int i = 0; i < 3; ++i)
42  {
43  if (_fe_problem.isTransient() && i < _ndisp)
44  _grad_disp_old[i] = &coupledGradientOld("displacements", i);
45  else
46  _grad_disp_old[i] = &_grad_zero;
47  }
48 }

Referenced by ComputeAxisymmetric1DFiniteStrain::initialSetup(), ComputeAxisymmetricRZIncrementalStrain::initialSetup(), ComputeRSphericalIncrementalStrain::initialSetup(), ComputeRSphericalFiniteStrain::initialSetup(), and ComputeAxisymmetric1DIncrementalStrain::initialSetup().

◆ initQpStatefulProperties()

void ComputeIncrementalStrainBase::initQpStatefulProperties ( )
overrideprotectedvirtualinherited

Reimplemented from ComputeStrainBase.

Reimplemented in ComputeCosseratIncrementalSmallStrain.

Definition at line 51 of file ComputeIncrementalStrainBase.C.

52 {
53  _mechanical_strain[_qp].zero();
54  _total_strain[_qp].zero();
55  _deformation_gradient[_qp].setToIdentity();
56 }

Referenced by ComputeCosseratIncrementalSmallStrain::initQpStatefulProperties().

◆ subtractEigenstrainIncrementFromStrain()

void ComputeIncrementalStrainBase::subtractEigenstrainIncrementFromStrain ( RankTwoTensor strain)
protectedinherited

Definition at line 59 of file ComputeIncrementalStrainBase.C.

60 {
61  for (unsigned int i = 0; i < _eigenstrains.size(); ++i)
62  {
63  strain -= (*_eigenstrains[i])[_qp];
64  strain += (*_eigenstrains_old[i])[_qp];
65  }
66 }

Referenced by ComputeIncrementalSmallStrain::computeProperties(), ComputeCosseratIncrementalSmallStrain::computeQpProperties(), and computeQpStrain().

◆ validParams()

InputParameters ComputeFiniteStrain::validParams ( )
static

Definition at line 27 of file ComputeFiniteStrain.C.

28 {
29  InputParameters params = ComputeIncrementalStrainBase::validParams();
30  params.addClassDescription(
31  "Compute a strain increment and rotation increment for finite strains.");
32  params.addParam<MooseEnum>("decomposition_method",
34  "Methods to calculate the strain and rotation increments");
35  return params;
36 }

Referenced by Compute1DFiniteStrain::validParams(), ComputeRSphericalFiniteStrain::validParams(), and Compute2DFiniteStrain::validParams().

Member Data Documentation

◆ _base_name

const std::string ComputeStrainBase::_base_name
protectedinherited

Definition at line 44 of file ComputeStrainBase.h.

Referenced by ComputeStrainBase::ComputeStrainBase().

◆ _current_elem_volume

const Real& ComputeStrainBase::_current_elem_volume
protectedinherited

◆ _decomposition_method

const DecompMethod ComputeFiniteStrain::_decomposition_method
private

Definition at line 46 of file ComputeFiniteStrain.h.

Referenced by computeQpIncrements().

◆ _deformation_gradient

MaterialProperty<RankTwoTensor>& ComputeIncrementalStrainBase::_deformation_gradient
protectedinherited

◆ _disp

std::vector<const VariableValue *> ComputeStrainBase::_disp
protectedinherited

◆ _eigenstrain_names

std::vector<MaterialPropertyName> ComputeStrainBase::_eigenstrain_names
protectedinherited

◆ _eigenstrains

std::vector<const MaterialProperty<RankTwoTensor> *> ComputeStrainBase::_eigenstrains
protectedinherited

◆ _eigenstrains_old

std::vector<const MaterialProperty<RankTwoTensor> *> ComputeIncrementalStrainBase::_eigenstrains_old
protectedinherited

◆ _Fhat

std::vector<RankTwoTensor> ComputeFiniteStrain::_Fhat
protected

◆ _global_strain

const MaterialProperty<RankTwoTensor>* ComputeStrainBase::_global_strain
protectedinherited

Definition at line 53 of file ComputeStrainBase.h.

Referenced by ComputeSmallStrain::computeProperties(), and computeQpStrain().

◆ _grad_disp

std::vector<const VariableGradient *> ComputeStrainBase::_grad_disp
protectedinherited

◆ _grad_disp_old

std::vector<const VariableGradient *> ComputeIncrementalStrainBase::_grad_disp_old
protectedinherited

◆ _mechanical_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBase::_mechanical_strain
protectedinherited

◆ _mechanical_strain_old

const MaterialProperty<RankTwoTensor>& ComputeIncrementalStrainBase::_mechanical_strain_old
protectedinherited

◆ _ndisp

unsigned int ComputeStrainBase::_ndisp
protectedinherited

◆ _rotation_increment

MaterialProperty<RankTwoTensor>& ComputeIncrementalStrainBase::_rotation_increment
protectedinherited

◆ _strain_increment

MaterialProperty<RankTwoTensor>& ComputeIncrementalStrainBase::_strain_increment
protectedinherited

◆ _strain_rate

MaterialProperty<RankTwoTensor>& ComputeIncrementalStrainBase::_strain_rate
protectedinherited

◆ _total_strain

MaterialProperty<RankTwoTensor>& ComputeStrainBase::_total_strain
protectedinherited

◆ _total_strain_old

const MaterialProperty<RankTwoTensor>& ComputeIncrementalStrainBase::_total_strain_old
protectedinherited

◆ _volumetric_locking_correction

const bool ComputeStrainBase::_volumetric_locking_correction
protectedinherited

The documentation for this class was generated from the following files:
ComputeFiniteStrain::_decomposition_method
const DecompMethod _decomposition_method
Definition: ComputeFiniteStrain.h:46
ComputeIncrementalStrainBase::_strain_rate
MaterialProperty< RankTwoTensor > & _strain_rate
Definition: ComputeIncrementalStrainBase.h:38
ComputeStrainBase::_current_elem_volume
const Real & _current_elem_volume
Definition: ComputeStrainBase.h:56
ComputeFiniteStrain::DecompMethod::EigenSolution
ComputeIncrementalStrainBase::subtractEigenstrainIncrementFromStrain
void subtractEigenstrainIncrementFromStrain(RankTwoTensor &strain)
Definition: ComputeIncrementalStrainBase.C:59
ComputeIncrementalStrainBase::ComputeIncrementalStrainBase
ComputeIncrementalStrainBase(const InputParameters &parameters)
Definition: ComputeIncrementalStrainBase.C:22
ComputeIncrementalStrainBase::_mechanical_strain_old
const MaterialProperty< RankTwoTensor > & _mechanical_strain_old
Definition: ComputeIncrementalStrainBase.h:44
ComputeIncrementalStrainBase::_grad_disp_old
std::vector< const VariableGradient * > _grad_disp_old
Definition: ComputeIncrementalStrainBase.h:36
ComputeFiniteStrain::_Fhat
std::vector< RankTwoTensor > _Fhat
Definition: ComputeFiniteStrain.h:37
ComputeIncrementalStrainBase::_rotation_increment
MaterialProperty< RankTwoTensor > & _rotation_increment
Definition: ComputeIncrementalStrainBase.h:40
ComputeIncrementalStrainBase::_eigenstrains_old
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains_old
Definition: ComputeIncrementalStrainBase.h:47
ComputeStrainBase::_ndisp
unsigned int _ndisp
Coupled displacement variables.
Definition: ComputeStrainBase.h:40
ComputeStrainBase::_mechanical_strain
MaterialProperty< RankTwoTensor > & _mechanical_strain
Definition: ComputeStrainBase.h:46
ComputeStrainBase::initialSetup
void initialSetup() override
Definition: ComputeStrainBase.C:75
ComputeIncrementalStrainBase::_total_strain_old
const MaterialProperty< RankTwoTensor > & _total_strain_old
Definition: ComputeIncrementalStrainBase.h:45
ComputeFiniteStrain::computeQpStrain
virtual void computeQpStrain()
Definition: ComputeFiniteStrain.C:105
ComputeIncrementalStrainBase::_deformation_gradient
MaterialProperty< RankTwoTensor > & _deformation_gradient
Definition: ComputeIncrementalStrainBase.h:42
ComputeStrainBase::_global_strain
const MaterialProperty< RankTwoTensor > * _global_strain
Definition: ComputeStrainBase.h:53
ComputeFiniteStrain::computeQpIncrements
virtual void computeQpIncrements(RankTwoTensor &e, RankTwoTensor &r)
Definition: ComputeFiniteStrain.C:137
ComputeStrainBase::_eigenstrains
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains
Definition: ComputeStrainBase.h:51
ComputeFiniteStrain::DecompMethod
DecompMethod
Definition: ComputeFiniteStrain.h:40
RankTwoTensorTempl< Real >
ComputeStrainBase::_total_strain
MaterialProperty< RankTwoTensor > & _total_strain
Definition: ComputeStrainBase.h:48
ComputeFiniteStrain::DecompMethod::TaylorExpansion
ComputeStrainBase::_volumetric_locking_correction
const bool _volumetric_locking_correction
Definition: ComputeStrainBase.h:55
ComputeIncrementalStrainBase::validParams
static InputParameters validParams()
Definition: ComputeIncrementalStrainBase.C:16
ComputeIncrementalStrainBase::_strain_increment
MaterialProperty< RankTwoTensor > & _strain_increment
Definition: ComputeIncrementalStrainBase.h:39
ComputeFiniteStrain::decompositionType
static MooseEnum decompositionType()
Definition: ComputeFiniteStrain.C:17
ComputeStrainBase::_grad_disp
std::vector< const VariableGradient * > _grad_disp
Definition: ComputeStrainBase.h:42