36 #ifndef _EXTENDEDNEWMARKBETA_HPP_ 37 #define _EXTENDEDNEWMARKBETA_HPP_ 39 #include <Eigen/Dense> 40 #include <Eigen/SparseCore> 65 ExtendedNewmarkBeta(std::shared_ptr<Mesh> &mesh,
double TimeStep,
double mtol=1E-12,
double ktol=1E-12,
double ftol=1E-12);
81 void SetAlgorithm(std::shared_ptr<Algorithm> &algorithm);
107 bool ComputeNewStep(std::shared_ptr<Mesh> &mesh,
unsigned int k=0);
139 void ComputeEffectiveForce(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff,
double factor=1.00,
unsigned int k=0);
168 Eigen::SparseMatrix<double>
M;
171 Eigen::SparseMatrix<double>
C;
174 Eigen::SparseMatrix<double>
K;
177 Eigen::SparseMatrix<double>
G;
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.