b2api
B2000++ API Reference Manual, VERSION 4.6
 
Loading...
Searching...
No Matches
b2000::solver::TimeIntegrator< T > Class Template Referenceabstract

Interface definition for performing time steps. More...

#include "b2time_integrator.H"

Inheritance diagram for b2000::solver::TimeIntegrator< T >:
Collaboration diagram for b2000::solver::TimeIntegrator< T >:

Public Member Functions

virtual void SetAttribute (const Case &b2case)=0
 
virtual void InitializeTimeIntegrator (const b2linalg::Matrix< T, b2linalg::Mrectangle > &initial_dof, const double dt)=0
 Initialize the time integrator with the model to act on the initial condition, the initial point in time and the case definition for the current computation.
 
virtual void InitializeFieldsForThisTimeStep (b2linalg::Matrix< T, b2linalg::Mrectangle > &dof, const b2linalg::Matrix< T, b2linalg::Mrectangle > &dU)=0
 Update velocities and (possibly) also accelerations after each time step. This method will also save the values from the previous time step ("current" t_n+1) as the new values for the following time step ("new" t_n). When this method is called the displacements dof[0] are already up to date, but velocities dof[1] and accelerations dof[2] are not.
 
virtual void UpdateFields (b2linalg::Matrix< T, b2linalg::Mrectangle > &dof, const b2linalg::Matrix< T, b2linalg::Mrectangle > &dU)=0
 Update velocities and (possibly) also accelerations after each global iteration.
 
virtual EffectiveSymmetricSystem GetEffectiveElementSystem (const b2linalg::Matrix< T, b2linalg::Mpacked > &k, const b2linalg::Matrix< T, b2linalg::Mpacked > &d, const b2linalg::Matrix< T, b2linalg::Mpacked > &m, const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &dof, const b2linalg::Index &index) const =0
 Calculates the element effective system based on the concrete time integrator.
 
virtual EffectiveUnsymmetricSystem GetEffectiveElementSystem (const b2linalg::Matrix< T, b2linalg::Mrectangle > &k, const b2linalg::Matrix< T, b2linalg::Mrectangle > &d, const b2linalg::Matrix< T, b2linalg::Mrectangle > &m, const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &dof, const b2linalg::Index &index) const =0
 Calculates the element effective system based on the concrete time integrator.
 
- Public Member Functions inherited from b2000::Object
virtual const std::string & get_object_name () const
 
virtual ~Object ()
 

Protected Attributes

b2linalg::Matrix< T, b2linalg::Mrectangle > V_ {}
 Displacements.
 
b2linalg::Matrix< T, b2linalg::Mrectangle > Z_ {}
 Velocities.
 
b2linalg::Matrix< T, b2linalg::Mrectangle > dot_Z_ {}
 Accelerations.
 
int current_stage_number_ {1}
 For multi-stage methods such as Runge-Kutta methods ...
 
int max_stage_number_ {1}
 ... (do not confuse this with b2000 stages)
 
double dt {1}
 Step size.
 

Additional Inherited Members

- Static Public Attributes inherited from b2000::Object
static ObjectType type
 

Detailed Description

template<typename T>
class b2000::solver::TimeIntegrator< T >

Interface definition for performing time steps.

Member Function Documentation

◆ GetEffectiveElementSystem() [1/2]

template<typename T >
virtual EffectiveSymmetricSystem b2000::solver::TimeIntegrator< T >::GetEffectiveElementSystem ( const b2linalg::Matrix< T, b2linalg::Mpacked > &  k,
const b2linalg::Matrix< T, b2linalg::Mpacked > &  d,
const b2linalg::Matrix< T, b2linalg::Mpacked > &  m,
const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &  dof,
const b2linalg::Index &  index 
) const
pure virtual

Calculates the element effective system based on the concrete time integrator.

Parameters
[in]kThe symmetric element stiffness matrix.
[in]dThe symmetric element damping matrix.
[in]mThe symmetric element mass matrix.
[in]dofCurrent global displacements, velocities and accelerations.
[in]indexLocal to global node mapping.
Returns
The effective element system with the symmetric matrix.

Implemented in b2000::solver::NewmarkIntegrator< T >.

◆ GetEffectiveElementSystem() [2/2]

template<typename T >
virtual EffectiveUnsymmetricSystem b2000::solver::TimeIntegrator< T >::GetEffectiveElementSystem ( const b2linalg::Matrix< T, b2linalg::Mrectangle > &  k,
const b2linalg::Matrix< T, b2linalg::Mrectangle > &  d,
const b2linalg::Matrix< T, b2linalg::Mrectangle > &  m,
const b2linalg::Matrix< T, b2linalg::Mrectangle_constref > &  dof,
const b2linalg::Index &  index 
) const
pure virtual

Calculates the element effective system based on the concrete time integrator.

Parameters
[in]kThe unsymmetric element stiffness matrix.
[in]dThe unsymmetric element damping matrix.
[in]mThe unsymmetric element mass matrix.
[in]dofCurrent global displacements, velocities and accelerations.
[in]indexLocal to global node mapping.
Returns
The effective element system with the unsymmetric matrix.

Implemented in b2000::solver::NewmarkIntegrator< T >.

◆ InitializeFieldsForThisTimeStep()

template<typename T >
virtual void b2000::solver::TimeIntegrator< T >::InitializeFieldsForThisTimeStep ( b2linalg::Matrix< T, b2linalg::Mrectangle > &  dof,
const b2linalg::Matrix< T, b2linalg::Mrectangle > &  dU 
)
pure virtual

Update velocities and (possibly) also accelerations after each time step. This method will also save the values from the previous time step ("current" t_n+1) as the new values for the following time step ("new" t_n). When this method is called the displacements dof[0] are already up to date, but velocities dof[1] and accelerations dof[2] are not.

Parameters
[in,out]dofFields containing displacements, velocities and accelerations.
[in]dUWorking array containing increments of current iteration.

Implemented in b2000::solver::NewmarkIntegrator< T >.

◆ InitializeTimeIntegrator()

template<typename T >
virtual void b2000::solver::TimeIntegrator< T >::InitializeTimeIntegrator ( const b2linalg::Matrix< T, b2linalg::Mrectangle > &  initial_dof,
const double  dt 
)
pure virtual

Initialize the time integrator with the model to act on the initial condition, the initial point in time and the case definition for the current computation.

Parameters
[in,out]initial_dofThe initial values of the DOF's.
[in]dtThe initial step size for the first time step.

Implemented in b2000::solver::NewmarkIntegrator< T >.

◆ SetAttribute()

template<typename T >
virtual void b2000::solver::TimeIntegrator< T >::SetAttribute ( const Case b2case)
pure virtual

Set the attributes of the time stepping methods according to the parameters set in the given case definition.

Parameters
[in]b2caseThe current case.

Implemented in b2000::solver::NewmarkIntegrator< T >.

◆ UpdateFields()

template<typename T >
virtual void b2000::solver::TimeIntegrator< T >::UpdateFields ( b2linalg::Matrix< T, b2linalg::Mrectangle > &  dof,
const b2linalg::Matrix< T, b2linalg::Mrectangle > &  dU 
)
pure virtual

Update velocities and (possibly) also accelerations after each global iteration.

Parameters
[in,out]dofFields containing displacements, velocities and accelerations.
[in]dUWorking array containing increments of current iteration.

Implemented in b2000::solver::NewmarkIntegrator< T >.

Member Data Documentation

◆ dot_Z_

template<typename T >
b2linalg::Matrix<T, b2linalg::Mrectangle> b2000::solver::TimeIntegrator< T >::dot_Z_ {}
protected

Accelerations.

For multi-stage methods dot_Z_ contains the velocities of all stages for each DOF, i.e. dot_Z_.size1() == total_num_dof, dot_Z_.size2() == number of stages.

For multi-step methods dot_Z_ contains the velocities of all required former time steps for each DOF, i.e. dot_Z_.size1() == total_num_dof, dot_Z_.size2() == k, for a k-step method.

◆ V_

template<typename T >
b2linalg::Matrix<T, b2linalg::Mrectangle> b2000::solver::TimeIntegrator< T >::V_ {}
protected

Displacements.

For multi-stage methods V_ contains the displacements only at time t_n for each DOF, i.e. V_.size1() == total_num_dof, V_.size2() == 1. Same holds for other one-step methods (e.g. Newmark).

For multi-step methods V_ contains the displacements of all required former time steps for each DOF, i.e. V_.size1() == total_num_dof, V_.size2() == k, for a k-step method.

◆ Z_

template<typename T >
b2linalg::Matrix<T, b2linalg::Mrectangle> b2000::solver::TimeIntegrator< T >::Z_ {}
protected

Velocities.

For multi-stage methods Z_ contains the velocities of all stages for each DOF, i.e. Z_.size1() == total_num_dof, Z_.size2() == number of stages.

For multi-step methods Z_ contains the velocities of all required former time steps for each DOF, i.e. Z_.size1() == total_num_dof, Z_.size2() == k, for a k-step method.


The documentation for this class was generated from the following file: