Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
lin3DTruss2.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 339-340.
25 // Prentice-Hall, 1996.
26 //
27 // Description:
30 //------------------------------------------------------------------------------
31 
32 #ifndef _LIN3DTRUSS2_HPP_
33 #define _LIN3DTRUSS2_HPP_
34 
35 #include <map>
36 #include <memory>
37 #include <string>
38 #include <Eigen/Dense>
39 
40 #include "Node.hpp"
41 #include "Load.hpp"
42 #include "Material.hpp"
43 #include "Element.hpp"
44 #include "Damping.hpp"
45 
53 class lin3DTruss2 : public Element{
54 
55  public:
62  lin3DTruss2(const std::vector<unsigned int> nodes, std::unique_ptr<Material> &material, const double area);
63 
65  ~lin3DTruss2();
66 
69  void CommitState();
70 
73  void ReverseState();
74 
77  void InitialState();
78 
81  void UpdateState();
82 
87  void SetDomain(std::map<unsigned int, std::shared_ptr<Node> > &nodes);
88 
92  void SetDamping(const std::shared_ptr<Damping> &damping);
93 
96  std::vector<unsigned int> GetTotalDegreeOfFreedom() const;
97 
101  Eigen::MatrixXd GetStrain() const;
102 
106  Eigen::MatrixXd GetStress() const;
107 
111  Eigen::MatrixXd GetStrainRate() const;
112 
118  Eigen::MatrixXd GetStrainAt(double x3, double x2) const;
119 
125  Eigen::MatrixXd GetStressAt(double x3, double x2) const;
126 
131  Eigen::VectorXd GetVTKResponse(std::string response) const;
132 
135  double ComputeEnergy();
136 
141  Eigen::MatrixXd ComputeMassMatrix();
142 
147  Eigen::MatrixXd ComputeStiffnessMatrix();
148 
153  Eigen::MatrixXd ComputeDampingMatrix();
154 
159  Eigen::MatrixXd ComputePMLMatrix();
160 
165  Eigen::VectorXd ComputeInternalForces();
166 
171  Eigen::VectorXd ComputeInternalDynamicForces();
172 
179  Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr<Load> &surface, unsigned int face);
180 
187  Eigen::VectorXd ComputeBodyForces(const std::shared_ptr<Load> &body, unsigned int k=0);
188 
195  Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr<Load> &drm, unsigned int k);
196 
197  private:
199  double Lo;
200 
202  double A;
203 
205  std::shared_ptr<Damping> theDamping;
206 
208  std::vector<std::shared_ptr<Node> > theNodes;
209 
211  std::unique_ptr<Material> theMaterial;
212 
214  double ComputeLength() const;
215 
219  Eigen::MatrixXd ComputeLocalAxes() const;
220 
224  Eigen::MatrixXd ComputeTransformationAxes() const;
225 
228  Eigen::VectorXd ComputeStrain() const;
229 
232  Eigen::VectorXd ComputeStrainRate() const;
233 
237  Eigen::MatrixXd ComputeInitialStiffnessMatrix() const;
238 };
239 
240 #endif
Eigen::VectorXd ComputeInternalDynamicForces()
Compute the elastic, inertial, and viscous forces acting on the element.
Eigen::MatrixXd ComputeLocalAxes() const
Compute/update the local axis-1 of the element.
double ComputeEnergy()
Computes the element energy for a given deformation.
Eigen::MatrixXd GetStrainRate() const
Gets the material/section (generalised) strain-rate.
This file contains the "Load" class declarations, which defines a load that can act in a node or an e...
Eigen::VectorXd GetVTKResponse(std::string response) const
Gets the element internal response in VTK format for Paraview display.
double Lo
Length of the element.
Definition: lin3DTruss2.hpp:199
Eigen::MatrixXd GetStress() const
Gets the material/section (generalised) stress.
Eigen::VectorXd ComputeSurfaceForces(const std::shared_ptr< Load > &surface, unsigned int face)
Compute the surface forces acting on the element.
std::vector< std::shared_ptr< Node > > theNodes
The Element&#39;s Nodes.
Definition: lin3DTruss2.hpp:208
This file contains the abstract "Material object" declarations, which computes the strain...
Eigen::MatrixXd GetStressAt(double x3, double x2) const
Gets the material stress in section at coordinate (x3,x2).
Eigen::VectorXd ComputeStrain() const
Update strain in the element.
Eigen::MatrixXd GetStrainAt(double x3, double x2) const
Gets the material strain in section at coordinate (x3,x2).
lin3DTruss2(const std::vector< unsigned int > nodes, std::unique_ptr< Material > &material, const double area)
Creates a lin3DTruss2 in a finite element Mesh.
std::vector< unsigned int > GetTotalDegreeOfFreedom() const
Gets the list of total-degree of freedom of this Element.
Eigen::VectorXd ComputeStrainRate() const
Update strain rate in the element.
std::shared_ptr< Damping > theDamping
The Damping model.
Definition: lin3DTruss2.hpp:205
Eigen::MatrixXd ComputeTransformationAxes() const
Compute/update the tranformation matrix from local to global axis.
double ComputeLength() const
Compute the current length of the element.
double A
Area of cross-section.
Definition: lin3DTruss2.hpp:202
Eigen::VectorXd ComputeBodyForces(const std::shared_ptr< Load > &body, unsigned int k=0)
Compute the body forces acting on the element.
Virtual class for creating an element in a mesh.
Definition: Element.hpp:51
Eigen::MatrixXd GetStrain() const
Gets the material/section (generalised) strain.
Eigen::VectorXd ComputeInternalForces()
Compute the internal (elastic) forces acting on the element.
Eigen::MatrixXd ComputePMLMatrix()
Compute the PML history matrix using gauss-integration.
Eigen::MatrixXd ComputeDampingMatrix()
Compute the damping matrix of the element using gauss-integration.
void CommitState()
Save the material states in the element.
void UpdateState()
Update the material states in the element.
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.
This file contains the "Damping object" declarations to update damping matrix to account for Caughey-...
std::unique_ptr< Material > theMaterial
The Element&#39;s material.
Definition: lin3DTruss2.hpp:211
Eigen::MatrixXd ComputeStiffnessMatrix()
Compute the stiffness matrix of the element using gauss-integration.
Eigen::VectorXd ComputeDomainReductionForces(const std::shared_ptr< Load > &drm, unsigned int k)
Compute the domain reduction forces acting on the element.
void ReverseState()
Reverse the material/section states to previous converged state in this element.
~lin3DTruss2()
Destroys this lin3DTruss2 object.
void SetDamping(const std::shared_ptr< Damping > &damping)
Sets the damping model.
Eigen::MatrixXd ComputeInitialStiffnessMatrix() const
Compute the initial stiffness matrix of the element.
void SetDomain(std::map< unsigned int, std::shared_ptr< Node > > &nodes)
Sets the finite element dependance among objects.
Eigen::MatrixXd ComputeMassMatrix()
Compute the lumped/consistent mass matrix of the element.
Class for creating a 3D linearized two-node truss element in a mesh.
Definition: lin3DTruss2.hpp:53
void InitialState()
Brings the material/section state to its initial state in this element.
This file contains the "Element" object declarations, which defines an element in a finite element me...