Go to the documentation of this file.
12 #include "ADMaterial.h"
13 #include "libmesh/dense_matrix.h"
15 #define usingComputeIncrementalShellStrainMembers \
16 usingMaterialMembers; \
17 using ADComputeIncrementalShellStrain<compute_stage>::_t_points; \
18 using ADComputeIncrementalShellStrain<compute_stage>::_2d_points; \
19 using ADComputeIncrementalShellStrain<compute_stage>::_nodes; \
20 using ADComputeIncrementalShellStrain<compute_stage>::_ge; \
21 using ADComputeIncrementalShellStrain<compute_stage>::_ge_old; \
22 using ADComputeIncrementalShellStrain<compute_stage>::_J_map; \
23 using ADComputeIncrementalShellStrain<compute_stage>::_J_map_old; \
24 using ADComputeIncrementalShellStrain<compute_stage>::_node_normal; \
25 using ADComputeIncrementalShellStrain<compute_stage>::_node_normal_old; \
26 using ADComputeIncrementalShellStrain<compute_stage>::_strain_vector; \
27 using ADComputeIncrementalShellStrain<compute_stage>::_soln_vector; \
28 using ADComputeIncrementalShellStrain<compute_stage>::_B; \
29 using ADComputeIncrementalShellStrain<compute_stage>::_strain_increment; \
30 using ADComputeIncrementalShellStrain<compute_stage>::_total_strain; \
31 using ADComputeIncrementalShellStrain<compute_stage>::_total_strain_old; \
32 using ADComputeIncrementalShellStrain<compute_stage>::_v1; \
33 using ADComputeIncrementalShellStrain<compute_stage>::_v2; \
34 using ADComputeIncrementalShellStrain<compute_stage>::_dxyz_dxi; \
35 using ADComputeIncrementalShellStrain<compute_stage>::_dxyz_deta; \
36 using ADComputeIncrementalShellStrain<compute_stage>::_dxyz_dzeta; \
37 using ADComputeIncrementalShellStrain<compute_stage>::_thickness; \
38 using ADComputeIncrementalShellStrain<compute_stage>::_dphidxi_map; \
39 using ADComputeIncrementalShellStrain<compute_stage>::_dphideta_map; \
40 using ADComputeIncrementalShellStrain<compute_stage>::_phi_map; \
41 using ADComputeIncrementalShellStrain<compute_stage>::_g1_a; \
42 using ADComputeIncrementalShellStrain<compute_stage>::_g1_c; \
43 using ADComputeIncrementalShellStrain<compute_stage>::_g2_b; \
44 using ADComputeIncrementalShellStrain<compute_stage>::_g2_d; \
45 using ADComputeIncrementalShellStrain<compute_stage>::_soln_disp_index; \
46 using ADComputeIncrementalShellStrain<compute_stage>::_soln_rot_index; \
47 using ADComputeIncrementalShellStrain<compute_stage>::_sol_old; \
48 using ADComputeIncrementalShellStrain<compute_stage>::computeBMatrix; \
49 using ADComputeIncrementalShellStrain<compute_stage>::computeSolnVector; \
50 using ADComputeIncrementalShellStrain<compute_stage>::updateGVectors; \
51 using ADComputeIncrementalShellStrain<compute_stage>::updatedxyz; \
52 using ADComputeIncrementalShellStrain<compute_stage>::computeGMatrix
55 template <ComputeStage>
70 template <ComputeStage compute_stage>
186 std::vector<ADRealVectorValue>
_v1;
189 std::vector<ADRealVectorValue>
_v2;
195 std::vector<const MaterialProperty<DenseMatrix<Real>> *>
_B_old;
201 std::vector<const MaterialProperty<RankTwoTensor> *>
_ge_old;
218 const NumericVector<Number> *
const &
_sol;
std::unique_ptr< QGauss > _t_qrule
Quadrature rule in the out of plane direction.
virtual void initQpStatefulProperties() override
RankTwoTensor _unrotated_total_strain
std::vector< ADMaterialProperty(RealVectorValue) * > _dxyz_dxi
Derivative of global x, y and z w.r.t isoparametric coordinate xi.
const bool _large_strain
Flag to compute large strains.
ADDenseVector _soln_vector
Vector that stores the incremental solution at all the 20 DOFs in the 4 noded element.
virtual void computeSolnVector()
Computes the 20x1 soln vector and its derivatives for each shell element.
std::vector< unsigned int > _disp_num
Variable numbers corresponding to the displacement variables.
declareADValidParams(ADComputeIncrementalShellStrain)
RankTwoTensorTempl< DualReal > DualRankTwoTensor
virtual void updatedxyz()
Updates covariant vectors at each qp for finite rotations.
std::vector< ADMaterialProperty(RankTwoTensor) * > _total_strain
Total strain increment in the covariant coordinate system.
std::vector< const MaterialProperty< RankTwoTensor > * > _ge_old
Old ge matrix for elasticity tensor conversion.
std::vector< std::vector< Real > > _phi_map
Shape function value.
std::vector< const MaterialProperty< RealVectorValue > * > _dxyz_dxi_old
Old derivative of global x, y and z w.r.t isoparametric coordinate xi.
std::vector< ADMaterialProperty(RealVectorValue) * > _dxyz_deta
Derivative of global x, y and z w.r.t isoparametric coordinate eta.
const NumericVector< Number > & _sol_old
std::vector< ADMaterialProperty(Real) * > _J_map
Material property containing jacobian of transformation.
std::vector< ADRealVectorValue > _v2
First tangential vectors at nodes.
std::vector< std::vector< Real > > _dphidxi_map
Derivatives of shape functions w.r.t isoparametric coordinates xi.
std::vector< std::vector< unsigned int > > _soln_disp_index
Indices of solution vector corresponding to displacement DOFs in 3 directions at the 4 nodes.
virtual void computeBMatrix()
Computes the B matrix that connects strains to nodal displacements and rotations.
std::vector< ADMaterialProperty(DenseMatrix< Real >) * > _B
B_matrix for small strain.
const MaterialProperty< RealVectorValue > & _node_normal_old
Material property storing the old normal to the element at the 4 nodes.
unsigned int _ndisp
Number of coupled displacement variables.
virtual void computeGMatrix()
Computes the transformation matrix from natural coordinates to local cartesian coordinates for elasti...
ADMaterialProperty(RealVectorValue) &_node_normal
Material property storing the normal to the element at the 4 nodes. Stored as a material property for...
const VariableValue & _thickness
Coupled variable for the shell thickness.
std::vector< ADMaterialProperty(RankTwoTensor) * > _ge
ge matrix for elasticity tensor conversion
unsigned int _nrot
Number of coupled rotational variables.
std::vector< ADMaterialProperty(RankTwoTensor) * > _strain_increment
Strain increment in the covariant coordinate system.
ADComputeIncrementalShellStrain(const InputParameters ¶meters)
std::vector< const MaterialProperty< RealVectorValue > * > _dxyz_deta_old
Old derivative of global x, y and z w.r.t isoparametric coordinate eta.
std::vector< MaterialProperty< RankTwoTensor > * > _rotation_matrix
Rotation matrix material property.
std::vector< const MaterialProperty< Real > * > _J_map_old
Old material property containing jacobian of transformation.
ADRealVectorValue _x2
simulation variables
NonlinearSystemBase & _nonlinear_sys
Reference to the nonlinear system object.
std::vector< Point > _t_points
Quadrature points in the out of plane direction in isoparametric coordinate system.
std::vector< ADMaterialProperty(RealVectorValue) * > _dxyz_dzeta
Derivative of global x, y and z w.r.t isoparametric coordinate zeta.
std::vector< const MaterialProperty< DenseMatrix< Real > > * > _B_old
Old B_matrix for small strain.
ADDenseVector _strain_vector
Vector that stores the strain in the the 2 axial and 3 shear directions.
std::vector< ADRealVectorValue > _v1
First tangential vectors at nodes.
std::vector< const MaterialProperty< RealVectorValue > * > _dxyz_dzeta_old
Old derivative of global x, y and z w.r.t isoparametric coordinate zeta.
std::vector< Point > _2d_points
Quadrature points in the in-plane direction in isoparametric coordinate system.
std::vector< MaterialProperty< RankTwoTensor > * > _total_global_strain
Total strain in global coordinate system.
std::vector< std::vector< Real > > _dphideta_map
Derivatives of shape functions w.r.t isoparametric coordinates eta.
std::vector< const Node * > _nodes
Vector storing pointers to the nodes of the shell element.
std::vector< unsigned int > _rot_num
Variable numbers corresponding to the rotational variables.
std::vector< std::vector< unsigned int > > _soln_rot_index
Indices of solution vector corresponding to rotation DOFs in 2 directions at the 4 nodes.
std::vector< const MaterialProperty< RankTwoTensor > * > _total_strain_old
Old total strain increment in the covariant coordinate system.
virtual void computeProperties() override
virtual void computeNodeNormal()
Computes the node normal at each node.
const NumericVector< Number > *const & _sol
RankTwoTensorTempl< Real > RankTwoTensor
virtual void updateGVectors()
Updates the vectors required for shear locking computation for finite rotations.