Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
QuasiStatic.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]
25 //
26 // Description:
30 //------------------------------------------------------------------------------
31 
32 #ifndef _QUASISTATIC_HPP_
33 #define _QUASISTATIC_HPP_
34 
35 #include <Eigen/Dense>
36 #include <Eigen/SparseCore>
37 
38 #include "Mesh.hpp"
39 #include "LoadCombo.hpp"
40 #include "Assembler.hpp"
41 #include "Algorithm.hpp"
42 #include "Integrator.hpp"
43 
51 class QuasiStatic : public Integrator{
52 
53  public:
61  QuasiStatic(std::shared_ptr<Mesh> &mesh, double mtol=1E-12, double ktol=1E-12, double ftol=1E-12);
62 
64  ~QuasiStatic();
65 
69  void Initialize(std::shared_ptr<Mesh> &mesh);
70 
73  void SetLoadCombination(std::shared_ptr<LoadCombo> &combo);
74 
77  void SetAlgorithm(std::shared_ptr<Algorithm> &algorithm);
78 
82  const Eigen::VectorXd& GetDisplacements();
83 
87  const Eigen::VectorXd& GetVelocities();
88 
92  const Eigen::VectorXd& GetAccelerations();
93 
97  const Eigen::VectorXd& GetPMLHistoryVector();
98 
103  bool ComputeNewStep(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
104 
111  Eigen::VectorXd ComputeReactionForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
112 
117  Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0);
118 
125  void ComputeSupportMotionVector(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0);
126 
134  void ComputeEffectiveForce(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0);
135 
141  void ComputeEffectiveStiffness(std::shared_ptr<Mesh> &mesh, Eigen::SparseMatrix<double> &Keff);
142 
143  protected:
145  double TimeStep;
146 
148  Eigen::VectorXd U;
149 
151  Eigen::VectorXd V;
152 
154  Eigen::VectorXd A;
155 
157  Eigen::VectorXd Ubar;
158 
160  Eigen::VectorXd Fbar;
161 
163  Eigen::SparseMatrix<double> K;
164 
166  std::weak_ptr<Algorithm> theAlgorithm;
167 
169  std::unique_ptr<Assembler> theAssembler;
170 };
171 
172 #endif
Eigen::VectorXd V
Total initial/previous velocity.
Definition: QuasiStatic.hpp:151
void SetAlgorithm(std::shared_ptr< Algorithm > &algorithm)
Sets the integrator for this algorithm.
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 QuasiStatic integrator.
const Eigen::VectorXd & GetVelocities()
Gets the velocity vector.
std::unique_ptr< Assembler > theAssembler
The finite element assembler.
Definition: QuasiStatic.hpp:169
Eigen::VectorXd Ubar
Total previous pml history values.
Definition: QuasiStatic.hpp:157
bool ComputeNewStep(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Computes a new time step.
This file contains the "Load Combination" class declarations, which defines how the loads are going t...
Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Gets the external force vector from previous analysis.
Class for defining an static integrator that neglects any damping and inertial forces present in the ...
Definition: QuasiStatic.hpp:51
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 > K
Model stiffness matrix.
Definition: QuasiStatic.hpp:163
double TimeStep
Integration time step.
Definition: QuasiStatic.hpp:145
void SetLoadCombination(std::shared_ptr< LoadCombo > &combo)
Set the load combination.
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 QuasiStatic integrator.
std::weak_ptr< Algorithm > theAlgorithm
The static solver algorithm.
Definition: QuasiStatic.hpp:166
This file contains the "Integrator object" declarations, and defines how the dynamic solution between...
~QuasiStatic()
Destroys this QuasiStatic object.
This file contains the "Mesh object" declarations, which stores nodes, materials, elements...
Eigen::VectorXd ComputeReactionForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)
Gets the reaction force ins this step.
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.
QuasiStatic(std::shared_ptr< Mesh > &mesh, double mtol=1E-12, double ktol=1E-12, double ftol=1E-12)
Creates a QuasiStatic object.
const Eigen::VectorXd & GetPMLHistoryVector()
Gets the PML history vector.
Eigen::VectorXd A
Total previous acceleration.
Definition: QuasiStatic.hpp:154
void Initialize(std::shared_ptr< Mesh > &mesh)
Initialize model matrices.
Eigen::VectorXd Fbar
The previous stage Force vector.
Definition: QuasiStatic.hpp:160
Eigen::VectorXd U
Total initial/previous displacements.
Definition: QuasiStatic.hpp:148
const Eigen::VectorXd & GetDisplacements()
Gets the displacement vector.