Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
ExtendedNewmarkBeta.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] Fathi, A., Poursartip, B., & Kallivokas, L. F. (2015). Time‐domain hybrid
25 // formulations for wave simulations in three‐dimensional PML‐truncated
26 // heterogeneous media. International Journal for Numerical Methods in
27 // Engineering, 101(3), 165-198.
28 //
29 // Description:
34 //------------------------------------------------------------------------------
35 
36 #ifndef _EXTENDEDNEWMARKBETA_HPP_
37 #define _EXTENDEDNEWMARKBETA_HPP_
38 
39 #include <Eigen/Dense>
40 #include <Eigen/SparseCore>
41 
42 #include "Mesh.hpp"
43 #include "LoadCombo.hpp"
44 #include "Assembler.hpp"
45 #include "Algorithm.hpp"
46 #include "Integrator.hpp"
47 
56 
57  public:
65  ExtendedNewmarkBeta(std::shared_ptr<Mesh> &mesh, double TimeStep, double mtol=1E-12, double ktol=1E-12, double ftol=1E-12);
66 
69 
73  void Initialize(std::shared_ptr<Mesh> &mesh);
74 
77  void SetLoadCombination(std::shared_ptr<LoadCombo> &combo);
78 
81  void SetAlgorithm(std::shared_ptr<Algorithm> &algorithm);
82 
86  const Eigen::VectorXd& GetDisplacements();
87 
91  const Eigen::VectorXd& GetVelocities();
92 
96  const Eigen::VectorXd& GetAccelerations();
97 
101  const Eigen::VectorXd& GetPMLHistoryVector();
102 
107  bool ComputeNewStep(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
108 
115  Eigen::VectorXd ComputeReactionForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
116 
121  Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
122 
130  void ComputeSupportMotionVector(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0);
131 
139  void ComputeEffectiveForce(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0);
140 
146  void ComputeEffectiveStiffness(std::shared_ptr<Mesh> &mesh, Eigen::SparseMatrix<double> &Keff);
147 
148  protected:
150  double dt;
151 
153  Eigen::VectorXd U;
154 
156  Eigen::VectorXd V;
157 
159  Eigen::VectorXd A;
160 
162  Eigen::VectorXd Ubar;
163 
165  Eigen::VectorXd Fbar;
166 
168  Eigen::SparseMatrix<double> M;
169 
171  Eigen::SparseMatrix<double> C;
172 
174  Eigen::SparseMatrix<double> K;
175 
177  Eigen::SparseMatrix<double> G;
178 
180  std::weak_ptr<Algorithm> theAlgorithm;
181 
183  std::unique_ptr<Assembler> theAssembler;
184 };
185 
186 #endif
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.
const Eigen::VectorXd & GetAccelerations()
Gets the acceleration vector.
void SetLoadCombination(std::shared_ptr< LoadCombo > &combo)
Set the load combination.
std::weak_ptr< Algorithm > theAlgorithm
The static solver algorithm.
Definition: ExtendedNewmarkBeta.hpp:180
Eigen::VectorXd V
Total initial/previous velocity.
Definition: ExtendedNewmarkBeta.hpp:156
This file contains the "Load Combination" class declarations, which defines how the loads are going t...
Eigen::VectorXd A
Total previous acceleration.
Definition: ExtendedNewmarkBeta.hpp:159
Eigen::VectorXd Ubar
Total previous pml history values.
Definition: ExtendedNewmarkBeta.hpp:162
void Initialize(std::shared_ptr< Mesh > &mesh)
Initialize model matrices.
Eigen::SparseMatrix< double > K
Model stiffness matrix.
Definition: ExtendedNewmarkBeta.hpp:174
std::unique_ptr< Assembler > theAssembler
The finite element assembler.
Definition: ExtendedNewmarkBeta.hpp:183
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
void SetAlgorithm(std::shared_ptr< Algorithm > &algorithm)
Sets the integrator for this algorithm.
Eigen::SparseMatrix< double > C
Model damping matrix.
Definition: ExtendedNewmarkBeta.hpp:171
void ComputeEffectiveStiffness(std::shared_ptr< Mesh > &mesh, Eigen::SparseMatrix< double > &Keff)
Gets the effective stiffness associated to the NewmarkBeta integrator.
const Eigen::VectorXd & GetDisplacements()
Gets the displacement vector.
Eigen::SparseMatrix< double > M
Model mass matrix.
Definition: ExtendedNewmarkBeta.hpp:168
This file contains the "Integrator object" declarations, and defines how the dynamic solution between...
double dt
Integration time step.
Definition: ExtendedNewmarkBeta.hpp:150
Class for integrating the third order differential equation using the implicit Newmark method...
Definition: ExtendedNewmarkBeta.hpp:55
Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Gets the external force vector from previous analysis.
Eigen::VectorXd U
Total initial/previous displacements.
Definition: ExtendedNewmarkBeta.hpp:153
Eigen::VectorXd ComputeReactionForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Gets the reaction force ins this step.
This file contains the "Mesh object" declarations, which stores nodes, materials, elements...
bool ComputeNewStep(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Computes a new time step.
ExtendedNewmarkBeta(std::shared_ptr< Mesh > &mesh, double TimeStep, double mtol=1E-12, double ktol=1E-12, double ftol=1E-12)
Creates a NewmarkBeta object.
Eigen::VectorXd Fbar
The previous stage Force vector.
Definition: ExtendedNewmarkBeta.hpp:165
~ExtendedNewmarkBeta()
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.
const Eigen::VectorXd & GetVelocities()
Gets the velocity vector.
Eigen::SparseMatrix< double > G
Model PML history matrix.
Definition: ExtendedNewmarkBeta.hpp:177
const Eigen::VectorXd & GetPMLHistoryVector()
Gets the PML history vector.