www.mooseframework.org
Public Member Functions | Static Public Member Functions | Protected Member Functions | Protected Attributes | Private Types | Private Attributes | List of all members
ADComputeRSphericalFiniteStrain< compute_stage > Class Template Reference

ADComputeRSphericalFiniteStrain defines a strain increment and a rotation increment for finite strains in 1D spherical symmetry geometries. More...

#include <ADComputeRSphericalFiniteStrain.h>

Inheritance diagram for ADComputeRSphericalFiniteStrain< compute_stage >:
[legend]

Public Member Functions

 ADComputeRSphericalFiniteStrain (const InputParameters &parameters)
 
virtual void initialSetup ()
 
virtual void computeProperties ()
 Computes the current and old deformation gradients with the assumptions for 1D spherical symmetry geometries: \( \epsilon_{\theta} = \epsilon_{\phi} = \frac{u_r}{r} \). More...
 

Static Public Member Functions

static InputParameters validParams ()
 
static MooseEnum decompositionType ()
 

Protected Member Functions

virtual void computeQpStrain ()
 
virtual void computeQpIncrements (ADRankTwoTensor &e, ADRankTwoTensor &r)
 
virtual void initQpStatefulProperties () override
 
void subtractEigenstrainIncrementFromStrain (ADRankTwoTensor &strain)
 
 ADMaterialProperty (RankTwoTensor) &_strain_rate
 
 ADMaterialProperty (RankTwoTensor) &_strain_increment
 
 ADMaterialProperty (RankTwoTensor) &_rotation_increment
 
const ADMaterialProperty (RankTwoTensor) *_global_strain
 
virtual void displacementIntegrityCheck ()
 

Protected Attributes

const VariableValue & _disp_old_0
 the old value of the first component of the displacements vector More...
 
 usingComputeFiniteStrainMembers
 
std::vector< ADRankTwoTensor > _Fhat
 
 usingComputeIncrementalStrainBaseMembers
 
std::vector< const VariableGradient * > _grad_disp_old
 
const MaterialProperty< RankTwoTensor > & _mechanical_strain_old
 
const MaterialProperty< RankTwoTensor > & _total_strain_old
 
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains_old
 
 usingComputeStrainBaseMembers
 
const unsigned int _ndisp
 Coupled displacement variables. More...
 
std::vector< const ADVariableValue * > _disp
 
std::vector< const ADVariableGradient * > _grad_disp
 
const std::string _base_name
 
std::vector< MaterialPropertyName > _eigenstrain_names
 
std::vector< const ADMaterialProperty(RankTwoTensor) * > _eigenstrains
 
const bool _volumetric_locking_correction
 
const Real & _current_elem_volume
 
 usingMaterialMembers
 

Private Types

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

Private Attributes

const DecompMethod _decomposition_method
 

Detailed Description

template<ComputeStage compute_stage>
class ADComputeRSphericalFiniteStrain< compute_stage >

ADComputeRSphericalFiniteStrain defines a strain increment and a rotation increment for finite strains in 1D spherical symmetry geometries.

The strains in the polar and azimuthal directions are functions of the radial displacement.

Definition at line 15 of file ADComputeRSphericalFiniteStrain.h.

Member Enumeration Documentation

◆ DecompMethod

template<ComputeStage compute_stage>
enum ADComputeFiniteStrain::DecompMethod
strongprivateinherited
Enumerator
TaylorExpansion 
EigenSolution 

Definition at line 47 of file ADComputeFiniteStrain.h.

48  {
49  TaylorExpansion,
50  EigenSolution
51  };

Constructor & Destructor Documentation

◆ ADComputeRSphericalFiniteStrain()

template<ComputeStage compute_stage>
ADComputeRSphericalFiniteStrain< compute_stage >::ADComputeRSphericalFiniteStrain ( const InputParameters &  parameters)

Definition at line 32 of file ADComputeRSphericalFiniteStrain.C.

35  _disp_old_0(coupledValueOld("displacements", 0))
36 {
37 }

Member Function Documentation

◆ ADMaterialProperty() [1/4]

template<ComputeStage compute_stage>
ADComputeIncrementalStrainBase< compute_stage >::ADMaterialProperty ( RankTwoTensor  ) &
protectedinherited

◆ ADMaterialProperty() [2/4]

template<ComputeStage compute_stage>
ADComputeIncrementalStrainBase< compute_stage >::ADMaterialProperty ( RankTwoTensor  ) &
protectedinherited

◆ ADMaterialProperty() [3/4]

template<ComputeStage compute_stage>
ADComputeIncrementalStrainBase< compute_stage >::ADMaterialProperty ( RankTwoTensor  ) &
protectedinherited

◆ ADMaterialProperty() [4/4]

template<ComputeStage compute_stage>
const ADComputeStrainBase< compute_stage >::ADMaterialProperty ( RankTwoTensor  )
protectedinherited

◆ computeProperties()

template<ComputeStage compute_stage>
void ADComputeRSphericalFiniteStrain< compute_stage >::computeProperties ( )
virtual

Computes the current and old deformation gradients with the assumptions for 1D spherical symmetry geometries: \( \epsilon_{\theta} = \epsilon_{\phi} = \frac{u_r}{r} \).

Definition at line 53 of file ADComputeRSphericalFiniteStrain.C.

54 {
55  for (_qp = 0; _qp < _qrule->n_points(); ++_qp)
56  {
57  // Deformation gradient calculation in cylindrical coordinates
58  ADRankTwoTensor A; // Deformation gradient
59  RankTwoTensor Fbar; // Old Deformation gradient
60 
61  // Step through calculating the current and old deformation gradients
62  // Only diagonal components are nonzero because this is a 1D material
63  // Note: x_disp is the radial displacement
64  A(0, 0) = (*_grad_disp[0])[_qp](0);
65  Fbar(0, 0) = (*_grad_disp_old[0])[_qp](0);
66 
67  // The polar and azimuthal strains are functions of radial displacement
68  if (!MooseUtils::relativeFuzzyEqual(_q_point[_qp](0), 0.0))
69  {
70  A(1, 1) = (*_disp[0])[_qp] / _q_point[_qp](0);
71  Fbar(1, 1) = _disp_old_0[_qp] / _q_point[_qp](0);
72  }
73 
74  // The polar and azimuthal strains are equivalent in this 1D problem
75  A(2, 2) = A(1, 1);
76  Fbar(2, 2) = Fbar(1, 1);
77 
78  // very nearly A = gradU - gradUold, adapted to cylindrical coords
79  A -= Fbar;
80 
81  // Fbar = ( I + gradUold)
82  Fbar.addIa(1.0);
83 
84  // Incremental deformation gradient _Fhat = I + A Fbar^-1
85  _Fhat[_qp] = A * Fbar.inverse();
86  _Fhat[_qp].addIa(1.0);
87 
89  }
90 }

◆ computeQpIncrements()

template<ComputeStage compute_stage>
void ADComputeFiniteStrain< compute_stage >::computeQpIncrements ( ADRankTwoTensor &  e,
ADRankTwoTensor &  r 
)
protectedvirtualinherited

Definition at line 128 of file ADComputeFiniteStrain.C.

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

◆ computeQpStrain()

template<ComputeStage compute_stage>
void ADComputeFiniteStrain< compute_stage >::computeQpStrain ( )
protectedvirtualinherited

Definition at line 95 of file ADComputeFiniteStrain.C.

96 {
97  ADRankTwoTensor total_strain_increment;
98 
99  // two ways to calculate these increments: TaylorExpansion(default) or EigenSolution
100  computeQpIncrements(total_strain_increment, _rotation_increment[_qp]);
101 
102  _strain_increment[_qp] = total_strain_increment;
103 
104  // Remove the eigenstrain increment
105  subtractEigenstrainIncrementFromStrain(_strain_increment[_qp]);
106 
107  if (_dt > 0)
108  _strain_rate[_qp] = _strain_increment[_qp] / _dt;
109  else
110  _strain_rate[_qp].zero();
111 
112  // Update strain in intermediate configuration
113  _mechanical_strain[_qp] = _mechanical_strain_old[_qp] + _strain_increment[_qp];
114  _total_strain[_qp] = _total_strain_old[_qp] + total_strain_increment;
115 
116  // Rotate strain to current configuration
117  _mechanical_strain[_qp] =
118  _rotation_increment[_qp] * _mechanical_strain[_qp] * _rotation_increment[_qp].transpose();
119  _total_strain[_qp] =
120  _rotation_increment[_qp] * _total_strain[_qp] * _rotation_increment[_qp].transpose();
121 
122  if (_global_strain)
123  _total_strain[_qp] += (*_global_strain)[_qp];
124 }

◆ decompositionType()

template<ComputeStage compute_stage>
MooseEnum ADComputeFiniteStrain< compute_stage >::decompositionType ( )
staticinherited

Definition at line 21 of file ADComputeFiniteStrain.C.

22 {
23  return MooseEnum("TaylorExpansion EigenSolution", "TaylorExpansion");
24 }

Referenced by ADComputeFiniteStrain< compute_stage >::validParams().

◆ displacementIntegrityCheck()

template<ComputeStage compute_stage>
void ADComputeStrainBase< compute_stage >::displacementIntegrityCheck ( )
protectedvirtualinherited

Reimplemented in ADCompute2DFiniteStrain< compute_stage >, ADCompute2DIncrementalStrain< compute_stage >, and ADCompute2DSmallStrain< compute_stage >.

Definition at line 93 of file ADComputeStrainBase.C.

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

◆ initialSetup()

template<ComputeStage compute_stage>
void ADComputeRSphericalFiniteStrain< compute_stage >::initialSetup ( )
virtual

Definition at line 41 of file ADComputeRSphericalFiniteStrain.C.

42 {
44 
45  const auto & subdomainIDs = _mesh.meshSubdomains();
46  for (auto subdomainID : subdomainIDs)
47  if (_fe_problem.getCoordSystem(subdomainID) != Moose::COORD_RSPHERICAL)
48  mooseError("The coordinate system must be set to RSPHERICAL for 1D R spherical simulations.");
49 }

◆ initQpStatefulProperties()

template<ComputeStage compute_stage>
void ADComputeIncrementalStrainBase< compute_stage >::initQpStatefulProperties ( )
overrideprotectedvirtualinherited

Reimplemented from ADComputeStrainBase< compute_stage >.

Definition at line 55 of file ADComputeIncrementalStrainBase.C.

56 {
57  _mechanical_strain[_qp].zero();
58  _total_strain[_qp].zero();
59 }

◆ subtractEigenstrainIncrementFromStrain()

template<ComputeStage compute_stage>
void ADComputeIncrementalStrainBase< compute_stage >::subtractEigenstrainIncrementFromStrain ( ADRankTwoTensor &  strain)
protectedinherited

Definition at line 63 of file ADComputeIncrementalStrainBase.C.

65 {
66  for (unsigned int i = 0; i < _eigenstrains.size(); ++i)
67  {
68  strain -= (*_eigenstrains[i])[_qp];
69  strain += (*_eigenstrains_old[i])[_qp];
70  }
71 }

◆ validParams()

template<ComputeStage compute_stage>
InputParameters ADComputeRSphericalFiniteStrain< compute_stage >::validParams ( )
static

Definition at line 23 of file ADComputeRSphericalFiniteStrain.C.

24 {
25  InputParameters params = ADComputeFiniteStrain<compute_stage>::validParams();
26  params.addClassDescription("Compute a strain increment and rotation increment for finite strains "
27  "in 1D spherical symmetry problems.");
28  return params;
29 }

Member Data Documentation

◆ _base_name

template<ComputeStage compute_stage>
const std::string ADComputeStrainBase< compute_stage >::_base_name
protectedinherited

◆ _current_elem_volume

template<ComputeStage compute_stage>
const Real& ADComputeStrainBase< compute_stage >::_current_elem_volume
protectedinherited

Definition at line 71 of file ADComputeStrainBase.h.

◆ _decomposition_method

template<ComputeStage compute_stage>
const DecompMethod ADComputeFiniteStrain< compute_stage >::_decomposition_method
privateinherited

Definition at line 53 of file ADComputeFiniteStrain.h.

◆ _disp

template<ComputeStage compute_stage>
std::vector<const ADVariableValue *> ADComputeStrainBase< compute_stage >::_disp
protectedinherited

Definition at line 57 of file ADComputeStrainBase.h.

◆ _disp_old_0

template<ComputeStage compute_stage>
const VariableValue& ADComputeRSphericalFiniteStrain< compute_stage >::_disp_old_0
protected

the old value of the first component of the displacements vector

Definition at line 41 of file ADComputeRSphericalFiniteStrain.h.

◆ _eigenstrain_names

template<ComputeStage compute_stage>
std::vector<MaterialPropertyName> ADComputeStrainBase< compute_stage >::_eigenstrain_names
protectedinherited

◆ _eigenstrains

template<ComputeStage compute_stage>
std::vector<const ADMaterialProperty(RankTwoTensor) *> ADComputeStrainBase< compute_stage >::_eigenstrains
protectedinherited

◆ _eigenstrains_old

template<ComputeStage compute_stage>
std::vector<const MaterialProperty<RankTwoTensor> *> ADComputeIncrementalStrainBase< compute_stage >::_eigenstrains_old
protectedinherited

◆ _Fhat

template<ComputeStage compute_stage>
std::vector<ADRankTwoTensor> ADComputeFiniteStrain< compute_stage >::_Fhat
protectedinherited

Definition at line 44 of file ADComputeFiniteStrain.h.

◆ _grad_disp

template<ComputeStage compute_stage>
std::vector<const ADVariableGradient *> ADComputeStrainBase< compute_stage >::_grad_disp
protectedinherited

Definition at line 58 of file ADComputeStrainBase.h.

◆ _grad_disp_old

template<ComputeStage compute_stage>
std::vector<const VariableGradient *> ADComputeIncrementalStrainBase< compute_stage >::_grad_disp_old
protectedinherited

Definition at line 49 of file ADComputeIncrementalStrainBase.h.

◆ _mechanical_strain_old

template<ComputeStage compute_stage>
const MaterialProperty<RankTwoTensor>& ADComputeIncrementalStrainBase< compute_stage >::_mechanical_strain_old
protectedinherited

Definition at line 55 of file ADComputeIncrementalStrainBase.h.

◆ _ndisp

template<ComputeStage compute_stage>
const unsigned int ADComputeStrainBase< compute_stage >::_ndisp
protectedinherited

Coupled displacement variables.

Definition at line 56 of file ADComputeStrainBase.h.

Referenced by ADComputeStrainBase< compute_stage >::ADComputeStrainBase().

◆ _total_strain_old

template<ComputeStage compute_stage>
const MaterialProperty<RankTwoTensor>& ADComputeIncrementalStrainBase< compute_stage >::_total_strain_old
protectedinherited

Definition at line 56 of file ADComputeIncrementalStrainBase.h.

◆ _volumetric_locking_correction

template<ComputeStage compute_stage>
const bool ADComputeStrainBase< compute_stage >::_volumetric_locking_correction
protectedinherited

◆ usingComputeFiniteStrainMembers

template<ComputeStage compute_stage>
ADComputeRSphericalFiniteStrain< compute_stage >::usingComputeFiniteStrainMembers
protected

Definition at line 43 of file ADComputeRSphericalFiniteStrain.h.

◆ usingComputeIncrementalStrainBaseMembers

template<ComputeStage compute_stage>
ADComputeFiniteStrain< compute_stage >::usingComputeIncrementalStrainBaseMembers
protectedinherited

Definition at line 56 of file ADComputeFiniteStrain.h.

◆ usingComputeStrainBaseMembers

template<ComputeStage compute_stage>
ADComputeIncrementalStrainBase< compute_stage >::usingComputeStrainBaseMembers
protectedinherited

Definition at line 60 of file ADComputeIncrementalStrainBase.h.

◆ usingMaterialMembers

template<ComputeStage compute_stage>
ADComputeStrainBase< compute_stage >::usingMaterialMembers
protectedinherited

Definition at line 73 of file ADComputeStrainBase.h.


The documentation for this class was generated from the following files:
ADComputeIncrementalStrainBase::_total_strain_old
const MaterialProperty< RankTwoTensor > & _total_strain_old
Definition: ADComputeIncrementalStrainBase.h:56
ADComputeFiniteStrain::DecompMethod::TaylorExpansion
ADComputeFiniteStrain::DecompMethod::EigenSolution
ADComputeIncrementalStrainBase::subtractEigenstrainIncrementFromStrain
void subtractEigenstrainIncrementFromStrain(ADRankTwoTensor &strain)
Definition: ADComputeIncrementalStrainBase.C:63
ADComputeFiniteStrain::_Fhat
std::vector< ADRankTwoTensor > _Fhat
Definition: ADComputeFiniteStrain.h:44
ADComputeFiniteStrain::_decomposition_method
const DecompMethod _decomposition_method
Definition: ADComputeFiniteStrain.h:53
ADComputeIncrementalStrainBase::_mechanical_strain_old
const MaterialProperty< RankTwoTensor > & _mechanical_strain_old
Definition: ADComputeIncrementalStrainBase.h:55
ADComputeFiniteStrain::validParams
static InputParameters validParams()
Definition: ADComputeFiniteStrain.C:28
ADComputeRSphericalFiniteStrain::_disp_old_0
const VariableValue & _disp_old_0
the old value of the first component of the displacements vector
Definition: ADComputeRSphericalFiniteStrain.h:41
ADComputeStrainBase::_eigenstrains
std::vector< const ADMaterialProperty(RankTwoTensor) * > _eigenstrains
Definition: ADComputeStrainBase.h:66
ADComputeFiniteStrain::computeQpStrain
virtual void computeQpStrain()
Definition: ADComputeFiniteStrain.C:95
ADComputeFiniteStrain
ADComputeFiniteStrain defines a strain increment and rotation increment, for finite strains.
Definition: ADComputeFiniteStrain.h:21
ADComputeStrainBase::_disp
std::vector< const ADVariableValue * > _disp
Definition: ADComputeStrainBase.h:57
ADComputeIncrementalStrainBase::_eigenstrains_old
std::vector< const MaterialProperty< RankTwoTensor > * > _eigenstrains_old
Definition: ADComputeIncrementalStrainBase.h:58
ADComputeIncrementalStrainBase::initialSetup
void initialSetup() override
Definition: ADComputeIncrementalStrainBase.C:41
ADComputeStrainBase::_grad_disp
std::vector< const ADVariableGradient * > _grad_disp
Definition: ADComputeStrainBase.h:58
RankTwoTensorTempl< Real >
ADComputeFiniteStrain::computeQpIncrements
virtual void computeQpIncrements(ADRankTwoTensor &e, ADRankTwoTensor &r)
Definition: ADComputeFiniteStrain.C:128
ADComputeIncrementalStrainBase::_grad_disp_old
std::vector< const VariableGradient * > _grad_disp_old
Definition: ADComputeIncrementalStrainBase.h:49
ADComputeStrainBase::_ndisp
const unsigned int _ndisp
Coupled displacement variables.
Definition: ADComputeStrainBase.h:56