37 #ifndef _NEWMARKBETA_HPP_ 38 #define _NEWMARKBETA_HPP_ 40 #include <Eigen/Dense> 41 #include <Eigen/SparseCore> 66 NewmarkBeta(std::shared_ptr<Mesh> &mesh,
double TimeStep,
double mtol=1E-12,
double ktol=1E-12,
double ftol=1E-12);
82 void SetAlgorithm(std::shared_ptr<Algorithm> &algorithm);
108 bool ComputeNewStep(std::shared_ptr<Mesh> &mesh,
unsigned int k=0);
140 void ComputeEffectiveForce(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff,
double factor=1.00,
unsigned int k=0);
169 Eigen::SparseMatrix<double>
M;
172 Eigen::SparseMatrix<double>
C;
175 Eigen::SparseMatrix<double>
K;
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