Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
lin3DFrame2.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] Finite Element Procedures, Bathe, K.J., Chapter 5: pages 400-408.
25 // Prentice-Hall, 1996.
26 // [2] Lars Andersen and Søren R.K. Nielsen, "Elastic Beams in Three Dimensions",
27 // Aalborg University, Department of Civil Engineering, ISSN 1901-7286,
28 // DCE Lecture Notes No. 23.
29 //
30 // Description:
33 //------------------------------------------------------------------------------
34 
35 #ifndef _LIN3DFRAME2_HPP_
36 #define _LIN3DFRAME2_HPP_
37 
38 #include <map>
39 #include <memory>
40 #include <string>
41 #include <Eigen/Dense>
42 
43 #include "Node.hpp"
44 #include "Load.hpp"
45 #include "Section.hpp"
46 #include "Element.hpp"
47 #include "Damping.hpp"
48 #include "QuadratureRule.hpp"
49 
57 class lin3DFrame2 : public Element{
58 
59  public:
68  lin3DFrame2(const std::vector<unsigned int> nodes, std::unique_ptr<Section> &section, bool formulation=false, const std::string quadrature="GAUSS", unsigned int nGauss=3);
69 
71  ~lin3DFrame2();
72 
75  void CommitState();
76 
79  void ReverseState();
80 
83  void InitialState();
84 
87  void UpdateState();
88 
93  void SetDomain(std::map<unsigned int, std::shared_ptr<Node> > &nodes);
94 
98  void SetDamping(const std::shared_ptr<Damping> &damping);
99 
102  std::vector<unsigned int> GetTotalDegreeOfFreedom() const;
103 
107  Eigen::MatrixXd GetStrain() const;
108 
112  Eigen::MatrixXd GetStress() const;
113 
117  Eigen::MatrixXd GetStrainRate() const;
118 
124  Eigen::MatrixXd GetStrainAt(double x3, double x2) const;
125 
131  Eigen::MatrixXd GetStressAt(double x3, double x2) const;
132 
137  Eigen::VectorXd GetVTKResponse(std::string response) const;
138 
141  double ComputeEnergy();
142 
147  Eigen::MatrixXd ComputeMassMatrix();
148 
153  Eigen::MatrixXd ComputeStiffnessMatrix();
154 
159  Eigen::MatrixXd ComputeDampingMatrix();
160 
165  Eigen::MatrixXd ComputePMLMatrix();
166 
171  Eigen::VectorXd ComputeInternalForces();
172 
177  Eigen::VectorXd ComputeInternalDynamicForces();
178 
185  Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr<Load> &surface, unsigned int face);
186 
193  Eigen::VectorXd ComputeBodyForces(const std::shared_ptr<Load> &body, unsigned int k=0);
194 
201  Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr<Load> &drm, unsigned int k);
202 
203  private:
205  double L;
206 
209 
211  double Phiy;
212 
214  double Phiz;
215 
217  std::shared_ptr<Damping> theDamping;
218 
220  std::vector<std::shared_ptr<Node> > theNodes;
221 
223  std::vector<std::unique_ptr<Section> > theSection;
224 
226  std::unique_ptr<QuadratureRule> QuadraturePoints;
227 
229  double ComputeLength() const;
230 
234  Eigen::MatrixXd ComputeLocalAxes() const;
235 
239  Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const;
240 
244  Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const;
245 
249  Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri) const;
250 
254  Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri) const;
255 
259  Eigen::MatrixXd ComputeInitialStiffnessMatrix() const;
260 };
261 
262 #endif
Eigen::MatrixXd ComputeStrainDisplacementMatrix(const double ri) const
Evaluates the strain-displacement matrix at a given Gauss point.
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
Eigen::MatrixXd ComputeShapeFunctionMatrix(const double ri) const
Evaluates the shape function matrix at a given Gauss point.
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
Class for creating a 3D linearized two-node frame element in a mesh.
Definition: lin3DFrame2.hpp:57
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element.
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
void InitialState()
Brings the material/section state to its initial state in this element.
void CommitState()
Save the material states in the element.
This file contains the abstract "Section object" declarations, which computes the strain...
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
double ComputeLength() const
Compute the current length of the element.
void UpdateState()
Update the material states in the element.
double Phiy
Section Shear Factor axis-2.
Definition: lin3DFrame2.hpp:211
Eigen::MatrixXd ComputeLocalAxes() const
Compute/update the local axis-1 of the element.
std::unique_ptr< QuadratureRule > QuadraturePoints
Coordinate of Gauss points.
Definition: lin3DFrame2.hpp:226
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
double Phiz
Section Shear Factor axis-3.
Definition: lin3DFrame2.hpp:214
~lin3DFrame2()
Destroys this lin3DFrame2 object.
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
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.
Eigen::VectorXd ComputeStrain(const Eigen::MatrixXd &Bij) const
Update strain in the element.
lin3DFrame2(const std::vector< unsigned int > nodes, std::unique_ptr< Section > &section, bool formulation=false, const std::string quadrature="GAUSS", unsigned int nGauss=3)
Creates a lin3DFrame2 in a finite element Mesh.
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
Eigen::MatrixXd ComputeInitialStiffnessMatrix() const
Compute the initial stiffness matrix of the element.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
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::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
This file contains the "Damping object" declarations to update damping matrix to account for Caughey-...
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
void ReverseState()
Reverse the material/section states to previous converged state in this element.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
bool Formulation
The element formulation.
Definition: lin3DFrame2.hpp:208
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: lin3DFrame2.hpp:217
Eigen::VectorXd ComputeStrainRate(const Eigen::MatrixXd &Bij) const
Update strain rate in the element.
std::vector< std::shared_ptr< Node > > theNodes
The Element&#39;s Nodes.
Definition: lin3DFrame2.hpp:220
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.
double L
Length of the element.
Definition: lin3DFrame2.hpp:205
double ComputeEnergy()
Computes the element energy for a given deformation.
This file contains the abstract "Quadrature Rule" declarations to integrate quantities in the iso-par...
std::vector< std::unique_ptr< Section > > theSection
The Element&#39;s Section.
Definition: lin3DFrame2.hpp:223
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
This file contains the "Element" object declarations, which defines an element in a finite element me...