Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
CompositeBathe.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 779-782.
25 // Prentice-Hall, 1996.
26 // [2] Conserving Energy and Momentum in Nonlinear Dynamics: A Simple Impicit
27 // Time Integration Scheme, Bathe, K.J., Computers and Structures, Vol(85),
28 // 437-445, 2007.
29 // [3] Insight into an implicit time integration scheme for structural dynamics,
30 // Klaus-Jurgen Bathe, Gunwoo Noh, Computers and Structures, Vol(98-99),
31 // 1-6, 2012.
32 //
33 // Description:
38 //------------------------------------------------------------------------------
39 
40 #ifndef _COMPOSITEBATHE_HPP_
41 #define _COMPOSITEBATHE_HPP_
42 
43 #include <Eigen/Dense>
44 #include <Eigen/SparseCore>
45 
46 #include "Mesh.hpp"
47 #include "LoadCombo.hpp"
48 #include "Assembler.hpp"
49 #include "Algorithm.hpp"
50 #include "Integrator.hpp"
51 
59 class CompositeBathe : public Integrator{
60 
61  public:
69  CompositeBathe(std::shared_ptr<Mesh> &mesh, double TimeStep, double mtol=1E-12, double ktol=1E-12, double ftol=1E-12);
70 
73 
77  void Initialize(std::shared_ptr<Mesh> &mesh);
78 
81  void SetLoadCombination(std::shared_ptr<LoadCombo> &combo);
82 
85  void SetAlgorithm(std::shared_ptr<Algorithm> &algorithm);
86 
90  const Eigen::VectorXd& GetDisplacements();
91 
95  const Eigen::VectorXd& GetVelocities();
96 
100  const Eigen::VectorXd& GetAccelerations();
101 
105  const Eigen::VectorXd& GetPMLHistoryVector();
106 
111  bool ComputeNewStep(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
112 
119  Eigen::VectorXd ComputeReactionForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
120 
125  Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
126 
134  void ComputeSupportMotionVector(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0);
135 
143  void ComputeEffectiveForce(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0);
144 
150  void ComputeEffectiveStiffness(std::shared_ptr<Mesh> &mesh, Eigen::SparseMatrix<double> &Keff);
151 
152  protected:
154  double dt;
155 
157  bool flag;
158 
160  Eigen::VectorXd U;
161 
163  Eigen::VectorXd V;
164 
166  Eigen::VectorXd A;
167 
169  Eigen::VectorXd Um;
170 
172  Eigen::VectorXd Vm;
173 
175  Eigen::VectorXd Ubar;
176 
178  Eigen::VectorXd Fbar;
179 
181  Eigen::SparseMatrix<double> M;
182 
184  Eigen::SparseMatrix<double> C;
185 
187  Eigen::SparseMatrix<double> K;
188 
190  std::weak_ptr<Algorithm> theAlgorithm;
191 
193  std::unique_ptr<Assembler> theAssembler;
194 };
195 
196 #endif
bool ComputeNewStep(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Computes a new time step.
void Initialize(std::shared_ptr< Mesh > &mesh)
Initialize model matrices.
void SetAlgorithm(std::shared_ptr< Algorithm > &algorithm)
Sets the integrator for this algorithm.
const Eigen::VectorXd & GetAccelerations()
Gets the acceleration vector.
void ComputeEffectiveStiffness(std::shared_ptr< Mesh > &mesh, Eigen::SparseMatrix< double > &Keff)
Gets the effective stiffness associated to the CompositeBathe integrator.
CompositeBathe(std::shared_ptr< Mesh > &mesh, double TimeStep, double mtol=1E-12, double ktol=1E-12, double ftol=1E-12)
Creates a CompositeBathe object.
void SetLoadCombination(std::shared_ptr< LoadCombo > &combo)
Set the load combination.
Class for integrating the equation of motion using an implicit second-order composite method...
Definition: CompositeBathe.hpp:59
This file contains the "Load Combination" class declarations, which defines how the loads are going t...
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.
Eigen::SparseMatrix< double > C
Model damping matrix.
Definition: CompositeBathe.hpp:184
Eigen::VectorXd Um
Total mid-point displacements.
Definition: CompositeBathe.hpp:169
Eigen::VectorXd U
Total current displacements.
Definition: CompositeBathe.hpp:160
This file contains the "Assembler object" declarations, this object computes the system mass...
Eigen::VectorXd ComputeReactionForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Gets the reaction force ins this step.
Class for defining how the static/dynamic analysis between two steps is going to be performed...
Definition: Integrator.hpp:53
Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Gets the external force vector from previous analysis.
const Eigen::VectorXd & GetPMLHistoryVector()
Gets the PML history vector.
This file contains the "Integrator object" declarations, and defines how the dynamic solution between...
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 CompositeBathe integrator.
const Eigen::VectorXd & GetVelocities()
Gets the velocity vector.
~CompositeBathe()
Destroys this CompositeBathe object.
This file contains the "Mesh object" declarations, which stores nodes, materials, elements...
const Eigen::VectorXd & GetDisplacements()
Gets the displacement vector.
Eigen::VectorXd V
Total current velocity.
Definition: CompositeBathe.hpp:163
std::weak_ptr< Algorithm > theAlgorithm
The static solver algorithm.
Definition: CompositeBathe.hpp:190
bool flag
Differentiate mid-point to end-point.
Definition: CompositeBathe.hpp:157
Eigen::VectorXd Vm
Total mid-point velocity.
Definition: CompositeBathe.hpp:172
Eigen::VectorXd Ubar
Total previous pml history values.
Definition: CompositeBathe.hpp:175
double dt
Integration time step.
Definition: CompositeBathe.hpp:154
std::unique_ptr< Assembler > theAssembler
The finite element assembler.
Definition: CompositeBathe.hpp:193
Eigen::SparseMatrix< double > K
Model stiffness matrix.
Definition: CompositeBathe.hpp:187
Eigen::VectorXd A
Total current acceleration.
Definition: CompositeBathe.hpp:166
Eigen::SparseMatrix< double > M
Model mass matrix.
Definition: CompositeBathe.hpp:181
Eigen::VectorXd Fbar
The previous stage Force vector.
Definition: CompositeBathe.hpp:178