KratosMultiphysics
KRATOS Multiphysics (Kratos) is a framework for building parallel, multi-disciplinary simulation software, aiming at modularity, extensibility, and high performance. Kratos is written in C++, and counts with an extensive Python interface.
Classes | List of all members
Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace > Class Template Reference

BDF integration scheme (for dynamic problems) More...

#include <residual_based_bdf_scheme.h>

Inheritance diagram for Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >:
Collaboration diagram for Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >:

Classes

struct  GeneralVectors
 

Public Member Functions

Life Cycle
 ResidualBasedBDFScheme (const std::size_t Order=2)
 Constructor. The BDF method. More...
 
 ResidualBasedBDFScheme (ResidualBasedBDFScheme &rOther)
 
BaseTypePointer Clone () override
 
 ~ResidualBasedBDFScheme () override
 
Operations
void Update (ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
 Performing the update of the solution. More...
 
void Predict (ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
 Performing the prediction of the solution. More...
 
void InitializeSolutionStep (ModelPart &rModelPart, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
 It initializes time step solution. Only for reasons if the time step solution is restarted. More...
 
int Check (const ModelPart &rModelPart) const override
 This function is designed to be called once to perform all the checks needed on the input provided. More...
 
void Clear () override
 Free memory allocated by this class. More...
 
Input and output
Parameters GetDefaultParameters () const override
 This method provides the defaults parameters to avoid conflicts between the different constructors. More...
 
std::string Info () const override
 Turn back information as a string. More...
 
void PrintInfo (std::ostream &rOStream) const override
 Print information about this object. More...
 
void PrintData (std::ostream &rOStream) const override
 Print object's data. More...
 
- Public Member Functions inherited from Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >
 ResidualBasedImplicitTimeScheme ()
 
 ResidualBasedImplicitTimeScheme (Parameters ThisParameters)
 Constructor. The implicit method method. More...
 
 ResidualBasedImplicitTimeScheme (ResidualBasedImplicitTimeScheme &rOther)
 
BaseType::Pointer Clone () override
 
 ~ResidualBasedImplicitTimeScheme () override
 
void CalculateSystemContributions (Element &rCurrentElement, LocalSystemMatrixType &rLHSContribution, LocalSystemVectorType &rRHSContribution, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
 This function is designed to be called in the builder and solver to introduce the selected time integration scheme. More...
 
void CalculateRHSContribution (Element &rCurrentElement, LocalSystemVectorType &rRHSContribution, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
 This function is designed to calculate just the RHS contribution. More...
 
void CalculateSystemContributions (Condition &rCurrentCondition, LocalSystemMatrixType &rLHSContribution, LocalSystemVectorType &rRHSContribution, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
 Functions totally analogous to the precedent but applied to the "condition" objects. More...
 
void CalculateRHSContribution (Condition &rCurrentCondition, LocalSystemVectorType &rRHSContribution, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
 Functions that calculates the RHS of a "condition" object. More...
 
void InitializeSolutionStep (ModelPart &rModelPart, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
 It initializes time step solution. Only for reasons if the time step solution is restarted. More...
 
int Check (const ModelPart &rModelPart) const override
 This function is designed to be called once to perform all the checks needed on the input provided. More...
 
Parameters GetDefaultParameters () const override
 This method provides the defaults parameters to avoid conflicts between the different constructors. More...
 
std::string Info () const override
 Turn back information as a string. More...
 
void PrintInfo (std::ostream &rOStream) const override
 Print information about this object. More...
 
void PrintData (std::ostream &rOStream) const override
 Print object's data. More...
 
 KRATOS_CLASS_POINTER_DEFINITION (ResidualBasedImplicitTimeScheme)
 
- Public Member Functions inherited from Kratos::Scheme< TSparseSpace, TDenseSpace >
 Scheme ()
 Default Constructor. More...
 
 Scheme (Parameters ThisParameters)
 Constructor with Parameters. More...
 
 Scheme (Scheme &rOther)
 
virtual ~Scheme ()
 
 KRATOS_CLASS_POINTER_DEFINITION (Scheme)
 Pointer definition of Scheme. More...
 
virtual ClassType::Pointer Create (Parameters ThisParameters) const
 Create method. More...
 
virtual void Initialize (ModelPart &rModelPart)
 This is the place to initialize the Scheme. More...
 
bool SchemeIsInitialized ()
 This method returns if the scheme is initialized. More...
 
void SetSchemeIsInitialized (bool SchemeIsInitializedFlag=true)
 This method sets if the elements have been initialized or not (true by default) More...
 
bool ElementsAreInitialized ()
 This method returns if the elements are initialized. More...
 
void SetElementsAreInitialized (bool ElementsAreInitializedFlag=true)
 This method sets if the elements have been initialized or not (true by default) More...
 
bool ConditionsAreInitialized ()
 This method returns if the conditions are initialized. More...
 
void SetConditionsAreInitialized (bool ConditionsAreInitializedFlag=true)
 This method sets if the conditions have been initialized or not (true by default) More...
 
virtual void InitializeElements (ModelPart &rModelPart)
 This is the place to initialize the elements. More...
 
virtual void InitializeConditions (ModelPart &rModelPart)
 This is the place to initialize the conditions. More...
 
virtual void FinalizeSolutionStep (ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
 Function called once at the end of a solution step, after convergence is reached if an iterative process is needed. More...
 
virtual void InitializeNonLinIteration (ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
 unction to be called when it is needed to initialize an iteration. It is designed to be called at the beginning of each non linear iteration More...
 
virtual void FinalizeNonLinIteration (ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
 Function to be called when it is needed to finalize an iteration. It is designed to be called at the end of each non linear iteration. More...
 
virtual void CalculateOutputData (ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
 Functions to be called to prepare the data needed for the output of results. More...
 
virtual void CleanOutputData ()
 Functions that cleans the results data. More...
 
virtual void Clean ()
 This function is intended to be called at the end of the solution step to clean up memory storage not needed after the end of the solution step. More...
 
virtual int Check (ModelPart &rModelPart)
 
virtual void CalculateLHSContribution (Element &rElement, LocalSystemMatrixType &LHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
 This function is designed to calculate just the LHS contribution. More...
 
virtual void CalculateLHSContribution (Condition &rCondition, LocalSystemMatrixType &LHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
 Functions totally analogous to the precedent but applied to the "condition" objects. More...
 
virtual void EquationId (const Element &rElement, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo)
 This method gets the eqaution id corresponding to the current element. More...
 
virtual void EquationId (const Condition &rCondition, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo)
 Functions totally analogous to the precedent but applied to the "condition" objects. More...
 
virtual void GetDofList (const Element &rElement, Element::DofsVectorType &rDofList, const ProcessInfo &rCurrentProcessInfo)
 Function that returns the list of Degrees of freedom to be assembled in the system for a Given element. More...
 
virtual void GetDofList (const Condition &rCondition, Element::DofsVectorType &rDofList, const ProcessInfo &rCurrentProcessInfo)
 Function that returns the list of Degrees of freedom to be assembled in the system for a Given condition. More...
 

Protected Member Functions

Protected Operations
void UpdateDerivatives (ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb)
 Performing the update of the derivatives. More...
 
virtual void UpdateFirstDerivative (NodesArrayType::iterator itNode)
 Updating first time derivative (velocity) More...
 
virtual void UpdateSecondDerivative (NodesArrayType::iterator itNode)
 Updating second time derivative (acceleration) More...
 
void AddDynamicsToLHS (LocalSystemMatrixType &rLHS_Contribution, LocalSystemMatrixType &rD, LocalSystemMatrixType &rM, const ProcessInfo &rCurrentProcessInfo) override
 It adds the dynamic LHS contribution of the elements. More...
 
template<class TObjectType >
void TemplateAddDynamicsToRHS (TObjectType &rObject, LocalSystemVectorType &rRHS_Contribution, LocalSystemMatrixType &rD, LocalSystemMatrixType &rM, const ProcessInfo &rCurrentProcessInfo)
 It adds the dynamic RHS contribution of the objects. More...
 
void AddDynamicsToRHS (Element &rElement, LocalSystemVectorType &rRHS_Contribution, LocalSystemMatrixType &rD, LocalSystemMatrixType &rM, const ProcessInfo &rCurrentProcessInfo) override
 It adds the dynamic RHS contribution of the elements. More...
 
void AddDynamicsToRHS (Condition &rCondition, LocalSystemVectorType &rRHS_Contribution, LocalSystemMatrixType &rD, LocalSystemMatrixType &rM, const ProcessInfo &rCurrentProcessInfo) override
 It adds the dynamic RHS contribution of the condition. More...
 
Protected Operations
- Protected Member Functions inherited from Kratos::Scheme< TSparseSpace, TDenseSpace >
virtual Parameters ValidateAndAssignParameters (Parameters ThisParameters, const Parameters DefaultParameters) const
 This method validate and assign default parameters. More...
 
virtual void AssignSettings (const Parameters ThisParameters)
 This method assigns settings to member variables. More...
 

Protected Attributes

Protected member Variables
const std::size_t mOrder
 
Vector mBDF
 The integration order. More...
 
GeneralVectors mVector
 The BDF coefficients. More...
 
Kratos::unique_ptr< TimeDiscretization::BDFmpBDFUtility
 The structure containing the derivatives. More...
 
- Protected Attributes inherited from Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >
GeneralMatrices mMatrix
 
- Protected Attributes inherited from Kratos::Scheme< TSparseSpace, TDenseSpace >
bool mSchemeIsInitialized
 
bool mElementsAreInitialized
 Flag to be used in controlling if the Scheme has been initialized or not. More...
 
bool mConditionsAreInitialized
 Flag taking in account if the elements were initialized correctly or not. More...
 

Type Definitions

typedef Scheme< TSparseSpace, TDenseSpace > BaseType
 
typedef BaseType::Pointer BaseTypePointer
 
typedef ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace > ImplicitBaseType
 
typedef ImplicitBaseType::TDataType TDataType
 
typedef ImplicitBaseType::DofsArrayType DofsArrayType
 
typedef Element::DofsVectorType DofsVectorType
 
typedef ImplicitBaseType::TSystemMatrixType TSystemMatrixType
 
typedef ImplicitBaseType::TSystemVectorType TSystemVectorType
 
typedef ImplicitBaseType::LocalSystemVectorType LocalSystemVectorType
 
typedef ImplicitBaseType::LocalSystemMatrixType LocalSystemMatrixType
 
typedef ModelPart::NodesContainerType NodesArrayType
 
static constexpr double ZeroTolerance = std::numeric_limits<double>::epsilon()
 Definition of epsilon. More...
 
 KRATOS_CLASS_POINTER_DEFINITION (ResidualBasedBDFScheme)
 

Additional Inherited Members

- Public Types inherited from Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >
typedef Scheme< TSparseSpace, TDenseSpace > BaseType
 Base class definition. More...
 
typedef BaseType::DofsArrayType DofsArrayType
 DoF array type definition. More...
 
typedef Element::DofsVectorType DofsVectorType
 DoF vector type definition. More...
 
typedef BaseType::TDataType TDataType
 Data type definition. More...
 
typedef BaseType::TSystemMatrixType TSystemMatrixType
 Matrix type definition. More...
 
typedef BaseType::TSystemVectorType TSystemVectorType
 Vector type definition. More...
 
typedef BaseType::LocalSystemVectorType LocalSystemVectorType
 Local system matrix type definition. More...
 
typedef BaseType::LocalSystemMatrixType LocalSystemMatrixType
 Local system vector type definition. More...
 
typedef ModelPart::NodesContainerType NodesArrayType
 Nodes containers definition. More...
 
typedef ModelPart::ElementsContainerType ElementsArrayType
 Elements containers definition. More...
 
typedef ModelPart::ConditionsContainerType ConditionsArrayType
 Conditions containers definition. More...
 
typedef std::size_t IndexType
 Index type definition. More...
 
- Public Types inherited from Kratos::Scheme< TSparseSpace, TDenseSpace >
using ClassType = Scheme< TSparseSpace, TDenseSpace >
 The definition of the current class. More...
 
using TDataType = typename TSparseSpace::DataType
 Data type definition. More...
 
using TSystemMatrixType = typename TSparseSpace::MatrixType
 Matrix type definition. More...
 
using TSystemVectorType = typename TSparseSpace::VectorType
 Vector type definition. More...
 
using LocalSystemMatrixType = typename TDenseSpace::MatrixType
 Local system matrix type definition. More...
 
using LocalSystemVectorType = typename TDenseSpace::VectorType
 Local system vector type definition. More...
 
using TDofType = Dof< double >
 DoF type definition. More...
 
using DofsArrayType = ModelPart::DofsArrayType
 DoF array type definition. More...
 
using ElementsArrayType = ModelPart::ElementsContainerType
 Elements containers definition. More...
 
using ConditionsArrayType = ModelPart::ConditionsContainerType
 Conditions containers definition. More...
 
- Static Public Member Functions inherited from Kratos::Scheme< TSparseSpace, TDenseSpace >
static std::string Name ()
 Returns the name of the class as used in the settings (snake_case format) More...
 

Detailed Description

template<class TSparseSpace, class TDenseSpace>
class Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >

BDF integration scheme (for dynamic problems)

The \( n \) order Backward Differentiation Formula (BDF) method is a two step \( n \) order accurate method. This scheme is designed to solve a system of the type:

\[ \mathbf{M} \frac{d^2(u_{n0})}{dt^2} + \mathbf{D} \frac{d(un0)}{dt} + \mathbf{K} u_{n0} = \mathbf{f}_{ext} \]

If we call:

Then we assume:

\[ \frac{d^2(u_{n0})}{dt^2} \|t_{n0} = \sum_i c_i \dot{u}_{ni} \]

\[ \frac{d(u_{n0})}{dt} \|t_{n0} = \sum_i c_i u_{n0} \]

with for order 2 (BDF2):

  1. \( c_0 = \frac{1.5}{dt} \)
  2. \( c_1 = \frac{-2.0}{dt} \)
  3. \( c_2 = \frac{0.5}{dt} \)

The LHS and RHS can be defined as:

\[ RHS = \mathbf{f}_{ext} - \mathbf{M} \frac{d(\dot{u}_{n0})}{dt} - \mathbf{D} \frac{d(u_{n0})}{dt} - \mathbf{K} u_{n0} \]

and

\[ LHS = \frac{d(-RHS)}{d(u_{n0})} = c_0^2 \mathbf{M} + c_0 \mathbf{D} + K \]

Note
This implies that elements are expected to be written in terms of a variable with two time derivatives Main reference
Todo:
Create a BibTeX file https://www.stack.nl/~dimitri/doxygen/manual/commands.html#cmdcite
Author
Vicente Mataix Ferrandiz

Member Typedef Documentation

◆ BaseType

template<class TSparseSpace , class TDenseSpace >
typedef Scheme<TSparseSpace,TDenseSpace> Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::BaseType

◆ BaseTypePointer

template<class TSparseSpace , class TDenseSpace >
typedef BaseType::Pointer Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::BaseTypePointer

◆ DofsArrayType

template<class TSparseSpace , class TDenseSpace >
typedef ImplicitBaseType::DofsArrayType Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::DofsArrayType

◆ DofsVectorType

template<class TSparseSpace , class TDenseSpace >
typedef Element::DofsVectorType Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::DofsVectorType

◆ ImplicitBaseType

template<class TSparseSpace , class TDenseSpace >
typedef ResidualBasedImplicitTimeScheme<TSparseSpace,TDenseSpace> Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::ImplicitBaseType

◆ LocalSystemMatrixType

template<class TSparseSpace , class TDenseSpace >
typedef ImplicitBaseType::LocalSystemMatrixType Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::LocalSystemMatrixType

◆ LocalSystemVectorType

template<class TSparseSpace , class TDenseSpace >
typedef ImplicitBaseType::LocalSystemVectorType Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::LocalSystemVectorType

◆ NodesArrayType

template<class TSparseSpace , class TDenseSpace >
typedef ModelPart::NodesContainerType Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::NodesArrayType

◆ TDataType

template<class TSparseSpace , class TDenseSpace >
typedef ImplicitBaseType::TDataType Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::TDataType

◆ TSystemMatrixType

template<class TSparseSpace , class TDenseSpace >
typedef ImplicitBaseType::TSystemMatrixType Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::TSystemMatrixType

◆ TSystemVectorType

template<class TSparseSpace , class TDenseSpace >
typedef ImplicitBaseType::TSystemVectorType Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::TSystemVectorType

Constructor & Destructor Documentation

◆ ResidualBasedBDFScheme() [1/2]

template<class TSparseSpace , class TDenseSpace >
Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::ResidualBasedBDFScheme ( const std::size_t  Order = 2)
inlineexplicit

Constructor. The BDF method.

Parameters
OrderThe integration order
Todo:
The ideal would be to use directly the dof or the variable itself to identify the type of variable and is derivatives

◆ ResidualBasedBDFScheme() [2/2]

template<class TSparseSpace , class TDenseSpace >
Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::ResidualBasedBDFScheme ( ResidualBasedBDFScheme< TSparseSpace, TDenseSpace > &  rOther)
inlineexplicit

Copy Constructor.

◆ ~ResidualBasedBDFScheme()

template<class TSparseSpace , class TDenseSpace >
Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::~ResidualBasedBDFScheme ( )
inlineoverride

Destructor.

Member Function Documentation

◆ AddDynamicsToLHS()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::AddDynamicsToLHS ( LocalSystemMatrixType rLHS_Contribution,
LocalSystemMatrixType rD,
LocalSystemMatrixType rM,
const ProcessInfo rCurrentProcessInfo 
)
inlineoverrideprotectedvirtual

It adds the dynamic LHS contribution of the elements.

\[ LHS = \frac{d(-RHS)}{d(u_{n0})} = c_0^2\mathbf{M} + c_0 \mathbf{D} + \mathbf{K} \]

Parameters
rLHS_ContributionThe dynamic contribution for the LHS
rDThe damping matrix
rMThe mass matrix
rCurrentProcessInfoThe current process info instance

Reimplemented from Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >.

Reimplemented in Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >.

◆ AddDynamicsToRHS() [1/2]

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::AddDynamicsToRHS ( Condition rCondition,
LocalSystemVectorType rRHS_Contribution,
LocalSystemMatrixType rD,
LocalSystemMatrixType rM,
const ProcessInfo rCurrentProcessInfo 
)
inlineoverrideprotectedvirtual

It adds the dynamic RHS contribution of the condition.

\[ RHS = f_{ext} - \ddot{u}_{n0} \mathbf{M} + \dot{u}_{n0} \mathbf{D} + u_{n0} \mathbf{K} \]

Parameters
rConditionThe condition to compute
RHS_ContributionThe dynamic contribution for the RHS
DThe damping matrix
MThe mass matrix
rCurrentProcessInfoThe current process info instance

Reimplemented from Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >.

Reimplemented in Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >.

◆ AddDynamicsToRHS() [2/2]

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::AddDynamicsToRHS ( Element rElement,
LocalSystemVectorType rRHS_Contribution,
LocalSystemMatrixType rD,
LocalSystemMatrixType rM,
const ProcessInfo rCurrentProcessInfo 
)
inlineoverrideprotectedvirtual

It adds the dynamic RHS contribution of the elements.

\[ \mathbf{b} - \mathbf{M} a - \mathbf{D} v \]

Parameters
rElementThe element to compute
RHS_ContributionThe dynamic contribution for the RHS
DThe damping matrix
MThe mass matrix
rCurrentProcessInfoThe current process info instance

Reimplemented from Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >.

Reimplemented in Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >.

◆ Check()

template<class TSparseSpace , class TDenseSpace >
int Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::Check ( const ModelPart rModelPart) const
inlineoverridevirtual

This function is designed to be called once to perform all the checks needed on the input provided.

Checks can be "expensive" as the function is designed to catch user's errors.

Parameters
rModelPartThe model of the problem to solve
Returns
Zero means all ok

Reimplemented from Kratos::Scheme< TSparseSpace, TDenseSpace >.

◆ Clear()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::Clear ( )
inlineoverridevirtual

Free memory allocated by this class.

Reimplemented from Kratos::Scheme< TSparseSpace, TDenseSpace >.

Reimplemented in Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >.

◆ Clone()

template<class TSparseSpace , class TDenseSpace >
BaseTypePointer Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::Clone ( )
inlineoverridevirtual

◆ GetDefaultParameters()

template<class TSparseSpace , class TDenseSpace >
Parameters Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::GetDefaultParameters ( ) const
inlineoverridevirtual

This method provides the defaults parameters to avoid conflicts between the different constructors.

Returns
The default parameters

Reimplemented from Kratos::Scheme< TSparseSpace, TDenseSpace >.

Reimplemented in Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >.

◆ Info()

template<class TSparseSpace , class TDenseSpace >
std::string Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::Info ( ) const
inlineoverridevirtual

◆ InitializeSolutionStep()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::InitializeSolutionStep ( ModelPart rModelPart,
TSystemMatrixType rA,
TSystemVectorType rDx,
TSystemVectorType rb 
)
inlineoverridevirtual

It initializes time step solution. Only for reasons if the time step solution is restarted.

Parameters
rModelPartThe model of the problem to solve
rALHS matrix
rDxIncremental update of primary variables
rbRHS Vector
Todo:
I cannot find the formula for the higher orders with variable time step. I tried to deduce by myself but the result was very unstable

Reimplemented from Kratos::Scheme< TSparseSpace, TDenseSpace >.

Reimplemented in Kratos::FluxCorrectedShallowWaterScheme< TSparseSpace, TDenseSpace >.

◆ KRATOS_CLASS_POINTER_DEFINITION()

template<class TSparseSpace , class TDenseSpace >
Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::KRATOS_CLASS_POINTER_DEFINITION ( ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >  )

◆ Predict()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::Predict ( ModelPart rModelPart,
DofsArrayType rDofSet,
TSystemMatrixType rA,
TSystemVectorType rDx,
TSystemVectorType rb 
)
inlineoverridevirtual

Performing the prediction of the solution.

It predicts the solution for the current step x = xold + vold * Dt

Parameters
rModelPartThe model of the problem to solve
rDofSetset of all primary variables
rALHS matrix
rDxIncremental update of primary variables
rbRHS Vector

Reimplemented from Kratos::Scheme< TSparseSpace, TDenseSpace >.

Reimplemented in Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >.

◆ PrintData()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::PrintData ( std::ostream &  rOStream) const
inlineoverridevirtual

Print object's data.

Reimplemented from Kratos::Scheme< TSparseSpace, TDenseSpace >.

◆ PrintInfo()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::PrintInfo ( std::ostream &  rOStream) const
inlineoverridevirtual

Print information about this object.

Reimplemented from Kratos::Scheme< TSparseSpace, TDenseSpace >.

◆ TemplateAddDynamicsToRHS()

template<class TSparseSpace , class TDenseSpace >
template<class TObjectType >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::TemplateAddDynamicsToRHS ( TObjectType &  rObject,
LocalSystemVectorType rRHS_Contribution,
LocalSystemMatrixType rD,
LocalSystemMatrixType rM,
const ProcessInfo rCurrentProcessInfo 
)
inlineprotected

It adds the dynamic RHS contribution of the objects.

\[ \mathbf{b} - \mathbf{M} a - \mathbf{D} v \]

Parameters
rObjectThe object to compute
rRHS_ContributionThe dynamic contribution for the RHS
rDThe damping matrix
rMThe mass matrix
rCurrentProcessInfoThe current process info instance

◆ Update()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::Update ( ModelPart rModelPart,
DofsArrayType rDofSet,
TSystemMatrixType rA,
TSystemVectorType rDx,
TSystemVectorType rb 
)
inlineoverridevirtual

Performing the update of the solution.

Incremental update within newton iteration. It updates the state variables at the end of the time step

\[ u_{n+1}^{k+1}= u_{n+1}^{k}+ \Delta u\]

Parameters
rModelPartThe model of the problem to solve
rDofSetSet of all primary variables
rALHS matrix
rDxincremental update of primary variables
rbRHS Vector

Reimplemented from Kratos::Scheme< TSparseSpace, TDenseSpace >.

Reimplemented in Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >.

◆ UpdateDerivatives()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::UpdateDerivatives ( ModelPart rModelPart,
DofsArrayType rDofSet,
TSystemMatrixType rA,
TSystemVectorType rDx,
TSystemVectorType rb 
)
inlineprotected

Performing the update of the derivatives.

Parameters
rModelPartThe model of the problem to solve
rDofSetSet of all primary variables
rALHS matrix
rDxincremental update of primary variables
rbRHS Vector

◆ UpdateFirstDerivative()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::UpdateFirstDerivative ( NodesArrayType::iterator  itNode)
inlineprotectedvirtual

◆ UpdateSecondDerivative()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::UpdateSecondDerivative ( NodesArrayType::iterator  itNode)
inlineprotectedvirtual

Member Data Documentation

◆ mBDF

template<class TSparseSpace , class TDenseSpace >
Vector Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::mBDF
protected

The integration order.

◆ mOrder

template<class TSparseSpace , class TDenseSpace >
const std::size_t Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::mOrder
protected

◆ mpBDFUtility

template<class TSparseSpace , class TDenseSpace >
Kratos::unique_ptr<TimeDiscretization::BDF> Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::mpBDFUtility
protected

The structure containing the derivatives.

◆ mVector

template<class TSparseSpace , class TDenseSpace >
GeneralVectors Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::mVector
protected

The BDF coefficients.

◆ ZeroTolerance

template<class TSparseSpace , class TDenseSpace >
constexpr double Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >::ZeroTolerance = std::numeric_limits<double>::epsilon()
staticconstexpr

Definition of epsilon.


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