Seismo-VLAB  1.3
An Open-Source Finite Element Software for Meso-Scale Simulations
Integrator Class Referenceabstract

Class for defining how the static/dynamic analysis between two steps is going to be performed. More...

#include <Integrator.hpp>

Inheritance diagram for Integrator:

Public Member Functions

 Integrator (const std::shared_ptr< Mesh > &mesh)
 Creates a Integrator object. More...
 
virtual ~Integrator ()=0
 Destroys this Integrator object. More...
 
virtual void Initialize (std::shared_ptr< Mesh > &mesh)=0
 Initialize model matrices. More...
 
virtual void SetLoadCombination (std::shared_ptr< LoadCombo > &combo)=0
 Set the load combination. More...
 
virtual void SetAlgorithm (std::shared_ptr< Algorithm > &algorithm)=0
 Sets the integrator for this algorithm. More...
 
virtual const Eigen::VectorXd & GetDisplacements ()=0
 Gets the displacement vector. More...
 
virtual const Eigen::VectorXd & GetVelocities ()=0
 Gets the velocity vector. More...
 
virtual const Eigen::VectorXd & GetAccelerations ()=0
 Gets the acceleration vector. More...
 
virtual const Eigen::VectorXd & GetPMLHistoryVector ()=0
 Gets the perfectly-matched layer history vector. More...
 
virtual bool ComputeNewStep (std::shared_ptr< Mesh > &mesh, unsigned int k=0)=0
 Computes a new time step. More...
 
virtual Eigen::VectorXd ComputeReactionForce (std::shared_ptr< Mesh > &mesh, unsigned int k=0)=0
 Gets the reaction force ins this step. More...
 
virtual Eigen::VectorXd ComputeProgressiveForce (std::shared_ptr< Mesh > &mesh, unsigned int k=0)=0
 Gets the external force vector from previous analysis. More...
 
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. More...
 
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. More...
 
virtual void ComputeEffectiveStiffness (std::shared_ptr< Mesh > &mesh, Eigen::SparseMatrix< double > &Keff)=0
 Gets the effective stiffness associated to the integrator. More...
 

Protected Attributes

Eigen::VectorXd SupportMotion
 Nodal support motion vector. More...
 
Eigen::SparseMatrix< double > Total2FreeMatrix
 Operator that enforced restrain/constraint. More...
 

Detailed Description

Class for defining how the static/dynamic analysis between two steps is going to be performed.

See also
QuasiStatic.hpp CentralDifference.hpp NewmarkBeta.hpp

Constructor & Destructor Documentation

◆ Integrator()

Integrator::Integrator ( const std::shared_ptr< Mesh > &  mesh)

Creates a Integrator object.

Parameters
meshPointer to the Mesh container to extract Node and Element.

◆ ~Integrator()

virtual Integrator::~Integrator ( )
pure virtual

Destroys this Integrator object.

Member Function Documentation

◆ ComputeEffectiveForce()

virtual void Integrator::ComputeEffectiveForce ( std::shared_ptr< Mesh > &  mesh,
Eigen::VectorXd &  Feff,
double  factor = 1.00,
unsigned int  k = 0 
)
pure virtual

Gets the effective force associated to the integrator.

Parameters
meshPointer to the Mesh object where Node and Element are stored.
FeffVector that stores the effective force.
factorThe incremental load factor.
kThe time step number to be solved.
See also
Assembler::ComputeInternalForceVector(), Assembler::ComputeExternalForceVector().

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ ComputeEffectiveStiffness()

virtual void Integrator::ComputeEffectiveStiffness ( std::shared_ptr< Mesh > &  mesh,
Eigen::SparseMatrix< double > &  Keff 
)
pure virtual

Gets the effective stiffness associated to the integrator.

Parameters
meshPointer to the Mesh object where Node and Element are stored.
KeffMatrix that stores the effective stiffness.
See also
Assembler::ComputeMassMatrix(), Assembler::ComputeStiffnessMatrix(), Assembler::ComputeDampingMatrix().

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ ComputeNewStep()

virtual bool Integrator::ComputeNewStep ( std::shared_ptr< Mesh > &  mesh,
unsigned int  k = 0 
)
pure virtual

Computes a new time step.

Parameters
meshThe finite element Mesh object.
kThe time step number to be solved.
Returns
Whether or not the Integrator was successfully applied.

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ ComputeProgressiveForce()

virtual Eigen::VectorXd Integrator::ComputeProgressiveForce ( std::shared_ptr< Mesh > &  mesh,
unsigned int  k = 0 
)
pure virtual

Gets the external force vector from previous analysis.

Parameters
meshPointer to the Mesh object where Node and Element are stored.
kThe time step number to be solved.
See also
Assembler::ComputeExternalForceVector().

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ ComputeReactionForce()

virtual Eigen::VectorXd Integrator::ComputeReactionForce ( std::shared_ptr< Mesh > &  mesh,
unsigned int  k = 0 
)
pure virtual

Gets the reaction force ins this step.

Parameters
meshPointer to the Mesh object where Node are stored.
kThe time step number to be solved.
Returns
Vector with the reaction forces and external forces.
Note
More details can be found at Reaction.
See also
Node::GetReaction(), Assembler::ComputeDynamicInternalForceVector().

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ ComputeSupportMotionVector()

virtual void Integrator::ComputeSupportMotionVector ( std::shared_ptr< Mesh > &  mesh,
Eigen::VectorXd &  Feff,
double  factor = 1.00,
unsigned int  k = 0 
)
pure virtual

Gets the incremental nodal support motion vector.

Parameters
meshPointer to the Mesh object where Node are stored.
FeffThe effective force vector to incorporate support motion forces.
factorThe incremental load factor.
kThe time step number to be solved.
See also
Node::GetSupportMotion(), Assembler::ComputeSupportMotionIncrement().

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ GetAccelerations()

virtual const Eigen::VectorXd& Integrator::GetAccelerations ( )
pure virtual

Gets the acceleration vector.

Returns
Vector with the acceleration states at current time step.

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ GetDisplacements()

virtual const Eigen::VectorXd& Integrator::GetDisplacements ( )
pure virtual

Gets the displacement vector.

Returns
Vector with the displacement states at current time step.

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ GetPMLHistoryVector()

virtual const Eigen::VectorXd& Integrator::GetPMLHistoryVector ( )
pure virtual

Gets the perfectly-matched layer history vector.

Returns
Vector with the PML history states at current time step.

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ GetVelocities()

virtual const Eigen::VectorXd& Integrator::GetVelocities ( )
pure virtual

Gets the velocity vector.

Returns
Vector with the velocity states at current time step.

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ Initialize()

virtual void Integrator::Initialize ( std::shared_ptr< Mesh > &  mesh)
pure virtual

Initialize model matrices.

Parameters
meshPointer to the Mesh container to extract Node and Element.
Note
This function computes matrices that are constant through the analysis.

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ SetAlgorithm()

virtual void Integrator::SetAlgorithm ( std::shared_ptr< Algorithm > &  algorithm)
pure virtual

Sets the integrator for this algorithm.

Parameters
algorithmPointer to the Algorithm to obtain the effective stiffness and force.

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

◆ SetLoadCombination()

virtual void Integrator::SetLoadCombination ( std::shared_ptr< LoadCombo > &  combo)
pure virtual

Set the load combination.

Parameters
comboPointer to the LoadCombo to be simulated.

Implemented in CompositeBathe, NewmarkBeta, ExtendedNewmarkBeta, CentralDifference, and QuasiStatic.

Member Data Documentation

◆ SupportMotion

Eigen::VectorXd Integrator::SupportMotion
protected

Nodal support motion vector.

◆ Total2FreeMatrix

Eigen::SparseMatrix<double> Integrator::Total2FreeMatrix
protected

Operator that enforced restrain/constraint.

See also
Mesh::GetTotalToFreeMatrix().