Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
lin3DShell4.hpp
Go to the documentation of this file.
1 //==============================================================================
2 //
3 // Seismo Virtual Laboratory
4 // Module for Serial and Parallel Analysis of seismic
5 // wave propagation and soil-structure interaction simulation
6 // Copyright (C) 2018-2021, The California Institute of Technology
7 // All Rights Reserved.
8 //
9 // Commercial use of this program without express permission of the California
10 // Institute of Technology, is strictly prohibited. See file "COPYRIGHT" in
11 // main directory for information on usage and redistribution, and for a
12 // DISCLAIMER OF ALL WARRANTIES.
13 //
14 //==============================================================================
15 //
16 // Written by:
17 // Danilo S. Kusanovic (dkusanov@caltech.edu)
18 // Elnaz E. Seylabi (elnaze@unr.edu)
19 //
20 // Supervised by:
21 // Domniki M. Asimaki (domniki@caltech.edu)
22 //
23 // References :
24 // [1] Adnan Ibrahimbegovic, "Quadrilateral finite elements for analysis of
25 // thick and thin plates", Computer Methods in Applied Mechanics and
26 // Engineering 110 (1993) 195-209.
27 // [2] Adnan Ibrahimbegovic, Robert Taylor, and Edward Wilson, "A Robust
28 // Quadrilateral Membrane Finite Element with Drilling degrees of freedom",
29 // International Journal for Numerical Methods in Engineering, vol. 30,
30 // 445-457 (1990)
31 //
32 // Description:
35 //------------------------------------------------------------------------------
36 
37 #ifndef _LIN3DSHELL4_HPP_
38 #define _LIN3DSHELL4_HPP_
39 
40 #include <map>
41 #include <memory>
42 #include <string>
43 #include <Eigen/Dense>
44 
45 #include "Node.hpp"
46 #include "Load.hpp"
47 #include "Section.hpp"
48 #include "Element.hpp"
49 #include "Damping.hpp"
50 #include "QuadratureRule.hpp"
51 
59 class lin3DShell4 : public Element{
60 
61  public:
69  lin3DShell4(const std::vector<unsigned int> nodes, std::unique_ptr<Section> &section, const std::string quadrature="GAUSS", const unsigned int nGauss=9);
70 
72  ~lin3DShell4();
73 
76  void CommitState();
77 
80  void ReverseState();
81 
84  void InitialState();
85 
88  void UpdateState();
89 
94  void SetDomain(std::map<unsigned int, std::shared_ptr<Node> > &nodes);
95 
99  void SetDamping(const std::shared_ptr<Damping> &damping);
100 
103  std::vector<unsigned int> GetTotalDegreeOfFreedom() const;
104 
108  Eigen::MatrixXd GetStrain() const;
109 
113  Eigen::MatrixXd GetStress() const;
114 
118  Eigen::MatrixXd GetStrainRate() const;
119 
125  Eigen::MatrixXd GetStrainAt(double x3, double x2) const;
126 
132  Eigen::MatrixXd GetStressAt(double x3, double x2) const;
133 
138  Eigen::VectorXd GetVTKResponse(std::string response) const;
139 
142  double ComputeEnergy();
143 
148  Eigen::MatrixXd ComputeMassMatrix();
149 
154  Eigen::MatrixXd ComputeStiffnessMatrix();
155 
160  Eigen::MatrixXd ComputeDampingMatrix();
161 
166  Eigen::MatrixXd ComputePMLMatrix();
167 
172  Eigen::VectorXd ComputeInternalForces();
173 
178  Eigen::VectorXd ComputeInternalDynamicForces();
179 
186  Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr<Load> &surface, unsigned int face);
187 
194  Eigen::VectorXd ComputeBodyForces(const std::shared_ptr<Load> &body, unsigned int k=0);
195 
202  Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr<Load> &drm, unsigned int k);
203 
204  private:
206  std::shared_ptr<Damping> theDamping;
207 
209  std::vector<std::shared_ptr<Node> > theNodes;
210 
212  std::vector<std::unique_ptr<Section> > theSection;
213 
215  std::unique_ptr<QuadratureRule> QuadraturePoints;
216 
218  double ComputeTotalMass();
219 
223  Eigen::MatrixXd ComputeLocalAxes() const;
224 
228  Eigen::MatrixXd ComputeRotation() const;
229 
232  Eigen::MatrixXd ComputeLocalCoordinates() const;
233 
239  Eigen::VectorXd ComputeStrain(Eigen::MatrixXd &BM12, Eigen::MatrixXd &xyloc, double ri, double si);
240 
246  Eigen::VectorXd ComputeStrainRate(Eigen::MatrixXd &BM12, Eigen::MatrixXd &xyloc, double ri, double si);
247 
249  Eigen::MatrixXd ComputeConstantTensionMatrix();
250 
255  void AssemblePlateMembraneEffects(Eigen::MatrixXd &A, const Eigen::MatrixXd &Am, const Eigen::MatrixXd &Ap);
256 
263  void ConstantTensionMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &BM12, Eigen::MatrixXd &Jij);
264 
271  void ComputeLoadShapeFunctionMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &Hij, Eigen::MatrixXd &Jij);
272 
279  void ComputePlateShapeFunctionMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &Hij, Eigen::MatrixXd &Jij);
280 
287  void ComputeMembraneShapeFunctionMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &Hij, Eigen::MatrixXd &Jij);
288 
295  void ComputePlateStrainDisplacementMatrix(const double ri, const double si, const Eigen::MatrixXd &xyloc, Eigen::MatrixXd &Bij, Eigen::MatrixXd &Jij);
296 
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);
306 
310  Eigen::MatrixXd ComputeInitialStiffnessMatrix();
311 };
312 
313 #endif
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 > &section, 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&#39;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&#39;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...