40 #ifndef _COMPOSITEBATHE_HPP_ 41 #define _COMPOSITEBATHE_HPP_ 43 #include <Eigen/Dense> 44 #include <Eigen/SparseCore> 69 CompositeBathe(std::shared_ptr<Mesh> &mesh,
double TimeStep,
double mtol=1E-12,
double ktol=1E-12,
double ftol=1E-12);
85 void SetAlgorithm(std::shared_ptr<Algorithm> &algorithm);
111 bool ComputeNewStep(std::shared_ptr<Mesh> &mesh,
unsigned int k=0);
143 void ComputeEffectiveForce(std::shared_ptr<Mesh> &mesh, Eigen::VectorXd &Feff,
double factor=1.00,
unsigned int k=0);
181 Eigen::SparseMatrix<double>
M;
184 Eigen::SparseMatrix<double>
C;
187 Eigen::SparseMatrix<double>
K;
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