Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
NewmarkBeta.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 6: pages 780-782.
25 // Prentice-Hall, 1996.
26 // [2] Consistent Finite-Element Response Sensitivity Analysis, J. P. Conte
27 // and P. K. Vijalapura and M. Meghella, Journal of Engineering Mechanics,
28 // volume 129 (12), 1380-1393, 2003.
29 //
30 // Description:
35 //------------------------------------------------------------------------------
36 
37 #ifndef _NEWMARKBETA_HPP_
38 #define _NEWMARKBETA_HPP_
39 
40 #include <Eigen/Dense>
41 #include <Eigen/SparseCore>
42 
43 #include "Mesh.hpp"
44 #include "LoadCombo.hpp"
45 #include "Assembler.hpp"
46 #include "Algorithm.hpp"
47 #include "Integrator.hpp"
48 
56 class NewmarkBeta : public Integrator{
57 
58  public:
66  NewmarkBeta(std::shared_ptr<Mesh> &mesh, double TimeStep, double mtol=1E-12, double ktol=1E-12, double ftol=1E-12);
67 
69  ~NewmarkBeta();
70 
74  void Initialize(std::shared_ptr<Mesh> &mesh);
75 
78  void SetLoadCombination(std::shared_ptr<LoadCombo> &combo);
79 
82  void SetAlgorithm(std::shared_ptr<Algorithm> &algorithm);
83 
87  const Eigen::VectorXd& GetDisplacements();
88 
92  const Eigen::VectorXd& GetVelocities();
93 
97  const Eigen::VectorXd& GetAccelerations();
98 
102  const Eigen::VectorXd& GetPMLHistoryVector();
103 
108  bool ComputeNewStep(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
109 
116  Eigen::VectorXd ComputeReactionForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
117 
122  Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
123 
131  void ComputeSupportMotionVector(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0);
132 
140  void ComputeEffectiveForce(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0);
141 
147  void ComputeEffectiveStiffness(std::shared_ptr<Mesh> &mesh, Eigen::SparseMatrix<double> &Keff);
148 
149  protected:
151  double dt;
152 
154  Eigen::VectorXd U;
155 
157  Eigen::VectorXd V;
158 
160  Eigen::VectorXd A;
161 
163  Eigen::VectorXd Ubar;
164 
166  Eigen::VectorXd Fbar;
167 
169  Eigen::SparseMatrix<double> M;
170 
172  Eigen::SparseMatrix<double> C;
173 
175  Eigen::SparseMatrix<double> K;
176 
178  std::weak_ptr<Algorithm> theAlgorithm;
179 
181  std::unique_ptr<Assembler> theAssembler;
182 };
183 
184 #endif
void Initialize(std::shared_ptr< Mesh > &mesh)
Initialize model matrices.
double dt
Integration time step.
Definition: NewmarkBeta.hpp:151
Eigen::SparseMatrix< double > C
Model damping matrix.
Definition: NewmarkBeta.hpp:172
~NewmarkBeta()
Destroys this NewmarkBeta object.
void ComputeEffectiveForce(std::shared_ptr< Mesh > &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0)
Gets the effective force associated to the NewmarkBeta integrator.
Eigen::VectorXd A
Total previous acceleration.
Definition: NewmarkBeta.hpp:160
Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Gets the external force vector from previous analysis.
Class for integrating the equation of motion using the implicit Newmark method.
Definition: NewmarkBeta.hpp:56
Eigen::VectorXd ComputeReactionForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Gets the reaction force ins this step.
void SetAlgorithm(std::shared_ptr< Algorithm > &algorithm)
Sets the integrator for this algorithm.
This file contains the "Load Combination" class declarations, which defines how the loads are going t...
bool ComputeNewStep(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Computes a new time step.
This file contains the "Assembler object" declarations, this object computes the system mass...
Class for defining how the static/dynamic analysis between two steps is going to be performed...
Definition: Integrator.hpp:53
Eigen::SparseMatrix< double > M
Model mass matrix.
Definition: NewmarkBeta.hpp:169
void SetLoadCombination(std::shared_ptr< LoadCombo > &combo)
Set the load combination.
const Eigen::VectorXd & GetAccelerations()
Gets the acceleration vector.
Eigen::VectorXd Ubar
Total previous pml history values.
Definition: NewmarkBeta.hpp:163
Eigen::VectorXd U
Total initial/previous displacements.
Definition: NewmarkBeta.hpp:154
const Eigen::VectorXd & GetVelocities()
Gets the velocity vector.
This file contains the "Integrator object" declarations, and defines how the dynamic solution between...
Eigen::VectorXd V
Total initial/previous velocity.
Definition: NewmarkBeta.hpp:157
void ComputeSupportMotionVector(std::shared_ptr< Mesh > &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0)
Gets the incremental nodal support motion vector.
std::weak_ptr< Algorithm > theAlgorithm
The static solver algorithm.
Definition: NewmarkBeta.hpp:178
This file contains the "Mesh object" declarations, which stores nodes, materials, elements...
NewmarkBeta(std::shared_ptr< Mesh > &mesh, double TimeStep, double mtol=1E-12, double ktol=1E-12, double ftol=1E-12)
Creates a NewmarkBeta object.
Eigen::SparseMatrix< double > K
Model stiffness matrix.
Definition: NewmarkBeta.hpp:175
std::unique_ptr< Assembler > theAssembler
The finite element assembler.
Definition: NewmarkBeta.hpp:181
void ComputeEffectiveStiffness(std::shared_ptr< Mesh > &mesh, Eigen::SparseMatrix< double > &Keff)
Gets the effective stiffness associated to the NewmarkBeta integrator.
const Eigen::VectorXd & GetPMLHistoryVector()
Gets the PML history vector.
const Eigen::VectorXd & GetDisplacements()
Gets the displacement vector.
Eigen::VectorXd Fbar
The previous stage Force vector.
Definition: NewmarkBeta.hpp:166