37 #ifndef _LIN3DSHELL4_HPP_ 38 #define _LIN3DSHELL4_HPP_ 43 #include <Eigen/Dense> 69 lin3DShell4(
const std::vector<unsigned int> nodes, std::unique_ptr<Section> §ion,
const std::string quadrature=
"GAUSS",
const unsigned int nGauss=9);
94 void SetDomain(std::map<
unsigned int, std::shared_ptr<Node> > &nodes);
99 void SetDamping(
const std::shared_ptr<Damping> &damping);
125 Eigen::MatrixXd
GetStrainAt(
double x3,
double x2)
const;
132 Eigen::MatrixXd
GetStressAt(
double x3,
double x2)
const;
194 Eigen::VectorXd
ComputeBodyForces(
const std::shared_ptr<Load> &body,
unsigned int k=0);
239 Eigen::VectorXd
ComputeStrain(Eigen::MatrixXd &BM12, Eigen::MatrixXd &xyloc,
double ri,
double si);
246 Eigen::VectorXd
ComputeStrainRate(Eigen::MatrixXd &BM12, Eigen::MatrixXd &xyloc,
double ri,
double si);
263 void ConstantTensionMatrix(
const double ri,
const double si,
const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &BM12, Eigen::MatrixXd &Jij);
305 void ComputeMembraneStrainDisplacementMatrix(
const double ri,
const double si,
const Eigen::MatrixXd &xyloc,
const Eigen::MatrixXd &Bij, Eigen::MatrixXd &BM12, Eigen::MatrixXd &Pij, Eigen::MatrixXd &Jij);
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
void ReverseState()
Reverse the material/section states to previous converged state in this element.
Class for creating a 3D linearized four-node shell element in a mesh.
Definition: lin3DShell4.hpp:59
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: lin3DShell4.hpp:215
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
lin3DShell4(const std::vector< unsigned int > nodes, std::unique_ptr< Section > §ion, const std::string quadrature="GAUSS", const unsigned int nGauss=9)
Creates a lin3DShell4 in a finite element Mesh.
Eigen::VectorXd ComputeStrain(Eigen::MatrixXd &BM12, Eigen::MatrixXd &xyloc, double ri, double si)
Update strain in the element.
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
void ComputePlateShapeFunctionMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &Hij, Eigen::MatrixXd &Jij)
Evaluates the shape function matrix for plate at a given Gauss point.
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
This file contains the abstract "Section object" declarations, which computes the strain...
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: lin3DShell4.hpp:206
~lin3DShell4()
Destroys this lin3DShell4 object.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
Eigen::VectorXd GetVTKResponse(std::string response) const
Gets the element internal response in VTK format for Paraview display.
void ComputeMembraneShapeFunctionMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &Hij, Eigen::MatrixXd &Jij)
Evaluates the shape function matrix for membrane at a given Gauss point.
void ComputeLoadShapeFunctionMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &Hij, Eigen::MatrixXd &Jij)
Evaluates the shape function matrix for plate at a given Gauss point.
Eigen::MatrixXd ComputeLocalCoordinates() const
Compute/update the node coordinates in axis-1-2-3 of the element.
Eigen::MatrixXd ComputeInitialStiffnessMatrix()
Compute the initial stiffness matrix of the element.
void ComputePlateStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &Bij, Eigen::MatrixXd &Jij)
Evaluates the strain-displacement matrix for a plate effect at a given Gauss point.
void InitialState()
Brings the material/section state to its initial state in this element.
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
void AssemblePlateMembraneEffects(Eigen::MatrixXd &A, const Eigen::MatrixXd &Am, const Eigen::MatrixXd &Ap)
Assembles the membrane and plate effects into a single matrix.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
double ComputeTotalMass()
Computes the total element mass.
double ComputeEnergy()
Computes the element energy for a given deformation.
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
void CommitState()
Save the material states in the element.
std::vector< std::unique_ptr< Section > > theSection
The Element's Sections.
Definition: lin3DShell4.hpp:212
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
This file contains the "Node object" declarations, which is stores the coordinates, state variables, degrees-of-freedom, and total-degree of freedom lists of a node in a finite element mesh.
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element.
This file contains the "Damping object" declarations to update damping matrix to account for Caughey-...
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
void ConstantTensionMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &BM12, Eigen::MatrixXd &Jij)
Evaluates the Constant Tension matrix and Jacobian for membrane effect at a given Gauss point...
Eigen::MatrixXd ComputeRotation() const
Compute/update the projection axis of the element.
void UpdateState()
Update the material states in the element.
Eigen::MatrixXd ComputeConstantTensionMatrix()
Correction matrix to Enforce Constant Tension for membrane effect (wilson).
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
Eigen::VectorXd ComputeStrainRate(Eigen::MatrixXd &BM12, Eigen::MatrixXd &xyloc, double ri, double si)
Update strain rate in the element.
void ComputeMembraneStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, const Eigen::MatrixXd &Bij, Eigen::MatrixXd &BM12, Eigen::MatrixXd &Pij, Eigen::MatrixXd &Jij)
Evaluates the strain-displacement matrix for membrane effect at a given Gauss point.
Eigen::MatrixXd ComputeLocalAxes() const
Compute/update the local axis-1-2-3 of the element.
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
std::vector< std::shared_ptr< Node > > theNodes
The Element's Nodes.
Definition: lin3DShell4.hpp:209
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
This file contains the "Element" object declarations, which defines an element in a finite element me...