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

ADStressDivergenceShell computes the stress divergence term for shell elements. More...

#include <ADStressDivergenceShell.h>

Inheritance diagram for ADStressDivergenceShell< compute_stage >:
[legend]

Public Member Functions

 ADStressDivergenceShell (const InputParameters &parameters)
 

Protected Member Functions

virtual ADReal computeQpResidual () override
 

Protected Attributes

const unsigned int _component
 
const bool _large_strain
 
std::vector< const ADMaterialProperty(RankTwoTensor) * > _stress
 
std::vector< const MaterialProperty< RankTwoTensor > * > _stress_old
 
std::vector< const ADMaterialProperty(DenseMatrix< Real >) * > _B_mat
 
std::vector< const ADMaterialProperty(DenseMatrix< Real >) * > _B_nl
 
std::vector< const ADMaterialProperty(Real) * > _J_map
 
std::unique_ptr< QGauss > _t_qrule
 Quadrature rule in the out of plane direction. More...
 
std::vector< Real > _t_weights
 Quadrature weights in the out of plane direction in isoparametric coordinate system. More...
 
std::vector< Real > _q_weights
 Qrule weights in isoparametric coordinate system. More...
 
unsigned int _qp_z
 qp index in out of plane direction More...
 
 usingKernelMembers
 

Detailed Description

template<ComputeStage compute_stage>
class ADStressDivergenceShell< compute_stage >

ADStressDivergenceShell computes the stress divergence term for shell elements.

Definition at line 19 of file ADStressDivergenceShell.h.

Constructor & Destructor Documentation

◆ ADStressDivergenceShell()

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

Definition at line 45 of file ADStressDivergenceShell.C.

46  : ADKernel<compute_stage>(parameters),
47  _component(getParam<unsigned int>("component")),
48  _large_strain(getParam<bool>("large_strain"))
49 {
50  _t_qrule = libmesh_make_unique<QGauss>(
51  1, Utility::string_to_enum<Order>(getParam<std::string>("through_thickness_order")));
52  _t_weights = _t_qrule->get_weights();
53 
54  _stress.resize(_t_weights.size());
55  _stress_old.resize(_t_weights.size());
56  _B_mat.resize(_t_weights.size());
57  if (_large_strain)
58  _B_nl.resize(_t_weights.size());
59  _J_map.resize(_t_weights.size());
60 
61  for (unsigned int i = 0; i < _t_weights.size(); ++i)
62  {
63  _stress[i] = &getADMaterialProperty<RankTwoTensor>("stress_t_points_" + std::to_string(i));
64  _stress_old[i] =
65  &getMaterialPropertyOldByName<RankTwoTensor>("stress_t_points_" + std::to_string(i));
66  _B_mat[i] = &getADMaterialProperty<DenseMatrix<Real>>("B_t_points_" + std::to_string(i));
67  if (_large_strain)
68  _B_nl[i] = &getADMaterialProperty<DenseMatrix<Real>>("B_nl_t_points_" + std::to_string(i));
69 
70  _J_map[i] = &getADMaterialProperty<Real>("J_mapping_t_points_" + std::to_string(i));
71  }
72 }

Member Function Documentation

◆ computeQpResidual()

template<ComputeStage compute_stage>
ADReal ADStressDivergenceShell< compute_stage >::computeQpResidual ( )
overrideprotectedvirtual

Definition at line 76 of file ADStressDivergenceShell.C.

77 {
78  _q_weights = _qrule->get_weights();
79  ADReal residual = 0.0;
80  ADReal residual1 = 0.0;
81  for (_qp_z = 0; _qp_z < _t_weights.size(); ++_qp_z)
82  {
83  residual1 = (*_stress[_qp_z])[_qp](0, 0) * (*_B_mat[_qp_z])[_qp](0, _i + _component * 4) +
84  (*_stress[_qp_z])[_qp](1, 1) * (*_B_mat[_qp_z])[_qp](1, _i + _component * 4) +
85  2.0 * (*_stress[_qp_z])[_qp](0, 1) * (*_B_mat[_qp_z])[_qp](2, _i + _component * 4) +
86  2.0 * (*_stress[_qp_z])[_qp](0, 2) * (*_B_mat[_qp_z])[_qp](3, _i + _component * 4) +
87  2.0 * (*_stress[_qp_z])[_qp](1, 2) * (*_B_mat[_qp_z])[_qp](4, _i + _component * 4);
88 
89  if (_large_strain)
90  residual1 +=
91  ((*_stress_old[_qp_z])[_qp](0, 0) * (*_B_nl[_qp_z])[_qp](0, _i + _component * 4) +
92  (*_stress_old[_qp_z])[_qp](1, 1) * (*_B_nl[_qp_z])[_qp](1, _i + _component * 4) +
93  2.0 * (*_stress_old[_qp_z])[_qp](0, 1) * (*_B_nl[_qp_z])[_qp](2, _i + _component * 4) +
94  2.0 * (*_stress_old[_qp_z])[_qp](0, 2) * (*_B_nl[_qp_z])[_qp](3, _i + _component * 4) +
95  2.0 * (*_stress_old[_qp_z])[_qp](1, 2) * (*_B_nl[_qp_z])[_qp](4, _i + _component * 4));
96 
97  residual += residual1 * (*_J_map[_qp_z])[_qp] * _q_weights[_qp] * _t_weights[_qp_z] /
98  (_ad_JxW[_qp] * _ad_coord[_qp]);
99  }
100  return residual;
101 }

Member Data Documentation

◆ _B_mat

template<ComputeStage compute_stage>
std::vector<const ADMaterialProperty(DenseMatrix<Real>) *> ADStressDivergenceShell< compute_stage >::_B_mat
protected

◆ _B_nl

template<ComputeStage compute_stage>
std::vector<const ADMaterialProperty(DenseMatrix<Real>) *> ADStressDivergenceShell< compute_stage >::_B_nl
protected

◆ _component

template<ComputeStage compute_stage>
const unsigned int ADStressDivergenceShell< compute_stage >::_component
protected

Definition at line 49 of file ADStressDivergenceShell.h.

◆ _J_map

template<ComputeStage compute_stage>
std::vector<const ADMaterialProperty(Real) *> ADStressDivergenceShell< compute_stage >::_J_map
protected

◆ _large_strain

template<ComputeStage compute_stage>
const bool ADStressDivergenceShell< compute_stage >::_large_strain
protected

◆ _q_weights

template<ComputeStage compute_stage>
std::vector<Real> ADStressDivergenceShell< compute_stage >::_q_weights
protected

Qrule weights in isoparametric coordinate system.

Definition at line 65 of file ADStressDivergenceShell.h.

◆ _qp_z

template<ComputeStage compute_stage>
unsigned int ADStressDivergenceShell< compute_stage >::_qp_z
protected

qp index in out of plane direction

Definition at line 68 of file ADStressDivergenceShell.h.

◆ _stress

template<ComputeStage compute_stage>
std::vector<const ADMaterialProperty(RankTwoTensor) *> ADStressDivergenceShell< compute_stage >::_stress
protected

◆ _stress_old

template<ComputeStage compute_stage>
std::vector<const MaterialProperty<RankTwoTensor> *> ADStressDivergenceShell< compute_stage >::_stress_old
protected

◆ _t_qrule

template<ComputeStage compute_stage>
std::unique_ptr<QGauss> ADStressDivergenceShell< compute_stage >::_t_qrule
protected

Quadrature rule in the out of plane direction.

Definition at line 59 of file ADStressDivergenceShell.h.

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

◆ _t_weights

template<ComputeStage compute_stage>
std::vector<Real> ADStressDivergenceShell< compute_stage >::_t_weights
protected

Quadrature weights in the out of plane direction in isoparametric coordinate system.

Definition at line 62 of file ADStressDivergenceShell.h.

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

◆ usingKernelMembers

template<ComputeStage compute_stage>
ADStressDivergenceShell< compute_stage >::usingKernelMembers
protected

Definition at line 70 of file ADStressDivergenceShell.h.


The documentation for this class was generated from the following files:
ADStressDivergenceShell::_stress_old
std::vector< const MaterialProperty< RankTwoTensor > * > _stress_old
Definition: ADStressDivergenceShell.h:53
ADStressDivergenceShell::_component
const unsigned int _component
Definition: ADStressDivergenceShell.h:49
ADStressDivergenceShell::_qp_z
unsigned int _qp_z
qp index in out of plane direction
Definition: ADStressDivergenceShell.h:68
ADStressDivergenceShell::_large_strain
const bool _large_strain
Definition: ADStressDivergenceShell.h:50
ADStressDivergenceShell::_J_map
std::vector< const ADMaterialProperty(Real) * > _J_map
Definition: ADStressDivergenceShell.h:56
ADStressDivergenceShell::_t_weights
std::vector< Real > _t_weights
Quadrature weights in the out of plane direction in isoparametric coordinate system.
Definition: ADStressDivergenceShell.h:62
ADStressDivergenceShell::_t_qrule
std::unique_ptr< QGauss > _t_qrule
Quadrature rule in the out of plane direction.
Definition: ADStressDivergenceShell.h:59
ADStressDivergenceShell::_B_nl
std::vector< const ADMaterialProperty(DenseMatrix< Real >) * > _B_nl
Definition: ADStressDivergenceShell.h:55
ADStressDivergenceShell::_q_weights
std::vector< Real > _q_weights
Qrule weights in isoparametric coordinate system.
Definition: ADStressDivergenceShell.h:65
ADStressDivergenceShell::_B_mat
std::vector< const ADMaterialProperty(DenseMatrix< Real >) * > _B_mat
Definition: ADStressDivergenceShell.h:54
ADStressDivergenceShell::_stress
std::vector< const ADMaterialProperty(RankTwoTensor) * > _stress
Definition: ADStressDivergenceShell.h:52