Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Integrator.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:
27 //This file contains the "Integrator object" declarations, and defines how the
28 //dynamic solution between two steps is going to be performed.
29 //------------------------------------------------------------------------------
30 
31 #ifndef _INTEGRATOR_HPP_
32 #define _INTEGRATOR_HPP_
33 
34 #include <memory>
35 #include <Eigen/Dense>
36 #include <Eigen/SparseCore>
37 
38 #include "Mesh.hpp"
39 #include "LoadCombo.hpp"
40 
41 class Algorithm;
42 
45 
53 class Integrator{
54 
55  public:
58  Integrator(const std::shared_ptr<Mesh> &mesh);
59 
61  virtual ~Integrator() = 0;
62 
66  virtual void Initialize(std::shared_ptr<Mesh> &mesh) = 0;
67 
70  virtual void SetLoadCombination(std::shared_ptr<LoadCombo> &combo) = 0;
71 
74  virtual void SetAlgorithm(std::shared_ptr<Algorithm> &algorithm) = 0;
75 
78  virtual const Eigen::VectorXd& GetDisplacements() = 0;
79 
82  virtual const Eigen::VectorXd& GetVelocities() = 0;
83 
86  virtual const Eigen::VectorXd& GetAccelerations() = 0;
87 
90  virtual const Eigen::VectorXd& GetPMLHistoryVector() = 0;
91 
96  virtual bool ComputeNewStep(std::shared_ptr<Mesh> &mesh, unsigned int k=0) = 0;
97 
104  virtual Eigen::VectorXd ComputeReactionForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0) = 0;
105 
110  virtual Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr<Mesh> &mesh, unsigned int k=0) = 0;
111 
118  virtual void ComputeSupportMotionVector(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0) = 0;
119 
126  virtual void ComputeEffectiveForce(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0) = 0;
127 
132  virtual void ComputeEffectiveStiffness(std::shared_ptr<Mesh> &mesh, Eigen::SparseMatrix<double> &Keff) = 0;
133 
134  protected:
136  Eigen::VectorXd SupportMotion;
137 
140  Eigen::SparseMatrix<double> Total2FreeMatrix;
141 
142  private:
143  //No private variables.
144 };
145 
146 #endif
Virtual class for defining how the solution between two steps is going to be carried out...
Definition: Algorithm.hpp:52
virtual const Eigen::VectorXd & GetPMLHistoryVector()=0
Gets the perfectly-matched layer history vector.
virtual const Eigen::VectorXd & GetVelocities()=0
Gets the velocity vector.
virtual Eigen::VectorXd ComputeProgressiveForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)=0
Gets the external force vector from previous analysis.
Eigen::SparseMatrix< double > Total2FreeMatrix
Operator that enforced restrain/constraint.
Definition: Integrator.hpp:140
virtual void ComputeSupportMotionVector(std::shared_ptr< Mesh > &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0)=0
Gets the incremental nodal support motion vector.
This file contains the "Load Combination" class declarations, which defines how the loads are going t...
virtual void SetAlgorithm(std::shared_ptr< Algorithm > &algorithm)=0
Sets the integrator for this algorithm.
virtual void ComputeEffectiveForce(std::shared_ptr< Mesh > &mesh, Eigen::VectorXd &Feff, double factor=1.00, unsigned int k=0)=0
Gets the effective force associated to the integrator.
Eigen::VectorXd SupportMotion
Nodal support motion vector.
Definition: Integrator.hpp:136
Class for defining how the static/dynamic analysis between two steps is going to be performed...
Definition: Integrator.hpp:53
virtual ~Integrator()=0
Destroys this Integrator object.
virtual Eigen::VectorXd ComputeReactionForce(std::shared_ptr< Mesh > &mesh, unsigned int k=0)=0
Gets the reaction force ins this step.
Integrator(const std::shared_ptr< Mesh > &mesh)
Creates a Integrator object.
virtual void ComputeEffectiveStiffness(std::shared_ptr< Mesh > &mesh, Eigen::SparseMatrix< double > &Keff)=0
Gets the effective stiffness associated to the integrator.
virtual const Eigen::VectorXd & GetDisplacements()=0
Gets the displacement vector.
This file contains the "Mesh object" declarations, which stores nodes, materials, elements...
virtual void Initialize(std::shared_ptr< Mesh > &mesh)=0
Initialize model matrices.
virtual bool ComputeNewStep(std::shared_ptr< Mesh > &mesh, unsigned int k=0)=0
Computes a new time step.
virtual const Eigen::VectorXd & GetAccelerations()=0
Gets the acceleration vector.
virtual void SetLoadCombination(std::shared_ptr< LoadCombo > &combo)=0
Set the load combination.