https://mooseframework.inl.gov
ADInertialForceShell.h
Go to the documentation of this file.
1 //* This file is part of the MOOSE framework
2 //* https://mooseframework.inl.gov
3 //*
4 //* All rights reserved, see COPYRIGHT for full restrictions
5 //* https://github.com/idaholab/moose/blob/master/COPYRIGHT
6 //*
7 //* Licensed under LGPL 2.1, please see LICENSE for details
8 //* https://www.gnu.org/licenses/lgpl-2.1.html
9 
10 #pragma once
11 
12 #include "ADTimeKernel.h"
13 #include "Material.h"
14 #include "RankTwoTensor.h"
15 #include "libmesh/dense_vector.h"
16 
17 namespace libMesh
18 {
19 class QGauss;
20 }
21 
23 {
24  std::array<ADRealVectorValue, 4> pos;
25  std::array<ADRealVectorValue, 4> rot;
26 };
27 
29 {
30 public:
32 
34 
35 protected:
36  virtual ADReal computeQpResidual() override { return 0.0; };
37 
38  virtual void computeResidual() override;
39  virtual void computeResidualsForJacobian() override;
41  const MooseArray<ADReal> & _ad_JxW);
42 
43 private:
45  unsigned int _nrot;
46 
48  unsigned int _ndisp;
49 
51  std::vector<unsigned int> _rot_num;
52 
54  std::vector<unsigned int> _disp_num;
55 
57  std::vector<unsigned int> _vel_num;
58 
60  std::vector<unsigned int> _accel_num;
61 
63  std::vector<unsigned int> _rot_vel_num;
64 
66  std::vector<unsigned int> _rot_accel_num;
67 
73 
75  const unsigned int _component;
76 
79 
82 
85 
88 
91 
94 
99  std::array<ADRealVectorValue, 4> _local_force, _local_moment;
100 
105  std::array<ADRealVectorValue, 4> _global_force, _global_moment;
106 
108  std::vector<std::vector<Real>> _dphidxi_map;
109 
111  std::vector<std::vector<Real>> _dphideta_map;
112 
114  std::vector<std::vector<Real>> _phi_map;
115 
117  std::vector<const Node *> _nodes;
118 
120  std::vector<Point> _2d_points;
121 
123  std::vector<ADReal> _2d_weights;
124 
126  std::vector<ADRealVectorValue> _v1;
127 
129  std::vector<ADRealVectorValue> _v2;
130 
133 
136 
138  std::vector<ADRealVectorValue> _node_normal;
139 
141  std::array<ADRealVectorValue, 2> _0g1_vectors;
142 
144  std::array<ADRealVectorValue, 2> _0g2_vectors;
145 
147  std::array<ADRealVectorValue, 2> _0g3_vectors;
148 
150  std::array<ADRealVectorValue, 2> _0g4_vectors;
151 
154 
157 
160 
163 
166 
168  const Real _alpha;
169 };
std::array< ADRealVectorValue, 2 > _0g2_vectors
Node 2 g vectors.
const Real _alpha
HHT time integration parameter.
std::vector< std::vector< Real > > _dphidxi_map
Derivatives of shape functions w.r.t isoparametric coordinates xi.
virtual void computeResidual() override
std::array< ADRealVectorValue, 2 > _0g3_vectors
Node 3 g vectors.
PosRotVectors _local_vel
Current shell nodal velocities in the local frame of reference.
PosRotVectors _old_vel
Old shell nodal velocities in the global frame of reference.
const MaterialProperty< Real > & _eta
Mass proportional Rayleigh damping parameter.
virtual ADReal computeQpResidual() override
PosRotVectors _local_accel
Current shell nodal accelerations in the local frame of reference.
const unsigned int _component
Direction along which residual is calculated.
std::array< ADRealVectorValue, 4 > _local_moment
std::array< ADRealVectorValue, 4 > _local_force
Forces and moments at the four nodes in the initial local configuration.
std::array< ADRealVectorValue, 4 > rot
std::vector< unsigned int > _vel_num
Variable numbers corresponding to velocity aux variables.
std::vector< unsigned int > _rot_vel_num
Variable numbers corresponding to rotational velocity aux variables.
PosRotVectors _local_old_vel
Old shell nodal velocities in the local frame of reference.
virtual void computeResidualsForJacobian() override
std::vector< unsigned int > _disp_num
Variable numbers corresponding to displacement variables.
The following methods are specializations for using the Parallel::packed_range_* routines for a vecto...
std::vector< ADRealVectorValue > _node_normal
Normal to the element at the 4 nodes.
std::vector< Point > _2d_points
Quadrature points in the in-plane direction in isoparametric coordinate system.
std::vector< ADReal > _2d_weights
Quadrature weights.
std::vector< std::vector< Real > > _dphideta_map
Derivatives of shape functions w.r.t isoparametric coordinates eta.
ADRealVectorValue _x3
Helper vector.
std::array< ADRealVectorValue, 2 > _0g1_vectors
Node 1 g vectors.
const ADMaterialProperty< RankTwoTensor > & _transformation_matrix
Rotation matrix material property.
std::vector< unsigned int > _accel_num
Variable numbers corresponding to acceleraion aux variables.
std::vector< const Node * > _nodes
Vector storing pointers to the nodes of the shell element.
std::vector< std::vector< Real > > _phi_map
Shape function value.
std::array< ADRealVectorValue, 4 > pos
const ADReal _thickness
Coupled variable for the shell thickness.
PosRotVectors _accel
Current shell nodal accelerations in the global frame of reference.
std::vector< unsigned int > _rot_num
Variable numbers corresponding to rotational variables.
std::array< ADRealVectorValue, 4 > _global_moment
PosRotVectors _vel
Current shell nodal velocities in the global frame of reference.
const MooseArray< ADReal > & _ad_JxW
std::vector< ADRealVectorValue > _v1
First tangential vectors at nodes.
virtual void computeShellInertialForces(const MooseArray< ADReal > &_ad_coord, const MooseArray< ADReal > &_ad_JxW)
DIE A HORRIBLE DEATH HERE typedef LIBMESH_DEFAULT_SCALAR_TYPE Real
std::array< ADRealVectorValue, 2 > _0g4_vectors
Node 4 g vectors.
std::vector< ADRealVectorValue > _v2
Second tangential vectors at nodes.
unsigned int _ndisp
Number of coupled displacement variables.
static InputParameters validParams()
unsigned int _nrot
Number of coupled rotational variables.
ADRankTwoTensor _original_local_config
Rotational transformation from global to initial shell local coordinate system.
const InputParameters & parameters() const
std::array< ADRealVectorValue, 4 > _global_force
Forces and moments at the four nodes in the global coordinate system.
const ADMaterialProperty< Real > & _J_map
Rotation matrix material property.
ADRealVectorValue _x2
Helper vector.
std::vector< unsigned int > _rot_accel_num
Variable numbers corresponding to rotational acceleration aux variables.
const MooseArray< ADReal > & _ad_coord
ADInertialForceShell(const InputParameters &parameters)
const MaterialProperty< Real > & _density
Shell material density.