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 | Protected Attributes | List of all members
Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace > Class Template Reference

Bossak integration scheme (for linear and nonlinear dynamic problems) for displacements. More...

#include <residual_based_bossak_displacement_scheme.hpp>

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

Classes

struct  BossakAlphaMethod
 Bossak Alpha parameters. More...
 
struct  GeneralVectors
 Velocities and accelerations used for integration. More...
 
struct  NewmarkMethod
 Newmark parameters used for integration. More...
 

Public Member Functions

Life Cycle
 ResidualBasedBossakDisplacementScheme (Parameters ThisParameters)
 Construct from a Parameters object. More...
 
 ResidualBasedBossakDisplacementScheme (const double Alpha=0.0)
 Constructor from a Bossak parameter. More...
 
 ResidualBasedBossakDisplacementScheme (const double Alpha, const double NewmarkBeta)
 Constructor. More...
 
 ResidualBasedBossakDisplacementScheme (ResidualBasedBossakDisplacementScheme &rOther)
 
BaseTypePointer Clone () override
 Clone method. More...
 
 ~ResidualBasedBossakDisplacementScheme () override
 
Input and output
std::string Info () const override
 Return 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 the instance'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 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 UpdateVelocity (array_1d< double, 3 > &rCurrentVelocity, const array_1d< double, 3 > &rDeltaDisplacement, const array_1d< double, 3 > &rPreviousVelocity, const array_1d< double, 3 > &rPreviousAcceleration)
 Update the first time derivative. More...
 
void UpdateAcceleration (array_1d< double, 3 > &rCurrentAcceleration, const array_1d< double, 3 > &rDeltaDisplacement, const array_1d< double, 3 > &rPreviousVelocity, const array_1d< double, 3 > &rPreviousAcceleration)
 Update the second time derivative. More...
 
void AddDynamicsToLHS (LocalSystemMatrixType &LHS_Contribution, LocalSystemMatrixType &D, LocalSystemMatrixType &M, const ProcessInfo &rCurrentProcessInfo) override
 Add dynamic left hand side contribution from Element s. More...
 
void AddDynamicsToRHS (Element &rElement, LocalSystemVectorType &RHS_Contribution, LocalSystemMatrixType &D, LocalSystemMatrixType &M, const ProcessInfo &rCurrentProcessInfo) override
 Add dynamic right hand side contribution from an Element. More...
 
void AddDynamicsToRHS (Condition &rCondition, LocalSystemVectorType &RHS_Contribution, LocalSystemMatrixType &D, LocalSystemMatrixType &M, const ProcessInfo &rCurrentProcessInfo) override
 Add dynamic right hand side contribution of a Condition. More...
 
void AssignSettings (const Parameters ThisParameters) override
 Assign member variables from Parameters. 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...
 

Protected Attributes

TSparseSpace::DofUpdaterPointerType mpDofUpdater = TSparseSpace::CreateDofUpdater()
 
BossakAlphaMethod mBossak
 Bossak Alpha parameters. More...
 
NewmarkMethod mNewmark
 Newmark Beta parameters. More...
 
GeneralVectors mVector
 Aggregate struct for velocities and accelerations. 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

using BaseType = Scheme< TSparseSpace, TDenseSpace >
 Base type for the scheme. More...
 
using ImplicitBaseType = ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >
 Implicit base type for the scheme. More...
 
using ClassType = ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >
 Class type for the scheme. More...
 
using TDataType = typename ImplicitBaseType::TDataType
 Data type used within the ImplicitBaseType. More...
 
using DofsArrayType = typename ImplicitBaseType::DofsArrayType
 Array type for degrees of freedom within ImplicitBaseType. More...
 
using DofsVectorType = typename Element::DofsVectorType
 Vector type for degrees of freedom within an Element. More...
 
using TSystemMatrixType = typename ImplicitBaseType::TSystemMatrixType
 Type for the system matrix within ImplicitBaseType. More...
 
using TSystemVectorType = typename ImplicitBaseType::TSystemVectorType
 Type for the system vector within ImplicitBaseType. More...
 
using LocalSystemVectorType = typename ImplicitBaseType::LocalSystemVectorType
 Type for local system vectors within ImplicitBaseType. More...
 
using LocalSystemMatrixType = typename ImplicitBaseType::LocalSystemMatrixType
 Type for local system matrices within ImplicitBaseType. More...
 
using NodeIterator = ModelPart::NodeIterator
 Iterator for nodes in a ModelPart. More...
 
using NodesArrayType = ModelPart::NodesContainerType
 Container type for nodes in a ModelPart. More...
 
using ElementsArrayType = ModelPart::ElementsContainerType
 Container type for elements in a ModelPart. More...
 
using ConditionsArrayType = ModelPart::ConditionsContainerType
 Container type for conditions in a ModelPart. More...
 
using BaseTypePointer = typename BaseType::Pointer
 Pointer type for the BaseType. More...
 
using ComponentType = double
 Component type as 'double'. More...
 
 KRATOS_CLASS_POINTER_DEFINITION (ResidualBasedBossakDisplacementScheme)
 

Operations

BaseType::Pointer Create (Parameters ThisParameters) const override
 Construct a dynamically allocated new scheme from a Parameters object. More...
 
void CalculateBossakCoefficients ()
 Recalculate the Newmark coefficients, taking the alpha parameters into account. More...
 
void Update (ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
 Update state variables within a newton iteration at the end of the time step. More...
 
void Predict (ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
 Apply the predictor. More...
 
void InitializeSolutionStep (ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b) override
 Prepare state variables for a new solution step. More...
 
int Check (const ModelPart &rModelPart) const override
 Check whether the scheme and the provided ModelPart are configured correctly. More...
 
void Clear () override
 Release dynamic memory allocated by this instance. More...
 
Parameters GetDefaultParameters () const override
 This function returns the default Parameters to help avoiding conflicts between different constructors. More...
 
static std::string Name ()
 Return the name of the class as used in the settings (snake_case). More...
 

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::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >

Bossak integration scheme (for linear and nonlinear dynamic problems) for displacements.

This is a dynamic implicit scheme based on the Bossak algorithm for displacements. The Bossak Alpha parameter ranges from 0 to -0.5, and introduces damping. Implementation according to: "An alpha modification of Newmark's method; W.L. Wood, M. Bossak, O.C. Zienkiewicz; Numerical Methods in Engineering; 1980"

Author
Josep Maria Carbonell
Vicente Mataix Ferrandiz
Andreas Winterstein (refactoring)

Member Typedef Documentation

◆ BaseType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::BaseType = Scheme<TSparseSpace, TDenseSpace>

Base type for the scheme.

◆ BaseTypePointer

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::BaseTypePointer = typename BaseType::Pointer

Pointer type for the BaseType.

◆ ClassType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::ClassType = ResidualBasedBossakDisplacementScheme<TSparseSpace, TDenseSpace>

Class type for the scheme.

◆ ComponentType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::ComponentType = double

Component type as 'double'.

◆ ConditionsArrayType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::ConditionsArrayType = ModelPart::ConditionsContainerType

Container type for conditions in a ModelPart.

◆ DofsArrayType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::DofsArrayType = typename ImplicitBaseType::DofsArrayType

Array type for degrees of freedom within ImplicitBaseType.

◆ DofsVectorType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::DofsVectorType = typename Element::DofsVectorType

Vector type for degrees of freedom within an Element.

◆ ElementsArrayType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::ElementsArrayType = ModelPart::ElementsContainerType

Container type for elements in a ModelPart.

◆ ImplicitBaseType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::ImplicitBaseType = ResidualBasedImplicitTimeScheme<TSparseSpace, TDenseSpace>

Implicit base type for the scheme.

◆ LocalSystemMatrixType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::LocalSystemMatrixType = typename ImplicitBaseType::LocalSystemMatrixType

Type for local system matrices within ImplicitBaseType.

◆ LocalSystemVectorType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::LocalSystemVectorType = typename ImplicitBaseType::LocalSystemVectorType

Type for local system vectors within ImplicitBaseType.

◆ NodeIterator

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::NodeIterator = ModelPart::NodeIterator

Iterator for nodes in a ModelPart.

◆ NodesArrayType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::NodesArrayType = ModelPart::NodesContainerType

Container type for nodes in a ModelPart.

◆ TDataType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::TDataType = typename ImplicitBaseType::TDataType

Data type used within the ImplicitBaseType.

◆ TSystemMatrixType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::TSystemMatrixType = typename ImplicitBaseType::TSystemMatrixType

Type for the system matrix within ImplicitBaseType.

◆ TSystemVectorType

template<class TSparseSpace , class TDenseSpace >
using Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::TSystemVectorType = typename ImplicitBaseType::TSystemVectorType

Type for the system vector within ImplicitBaseType.

Constructor & Destructor Documentation

◆ ResidualBasedBossakDisplacementScheme() [1/4]

template<class TSparseSpace , class TDenseSpace >
Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::ResidualBasedBossakDisplacementScheme ( Parameters  ThisParameters)
inlineexplicit

Construct from a Parameters object.

Parameters
ThisParametersThe parameters containing the configuration.

◆ ResidualBasedBossakDisplacementScheme() [2/4]

template<class TSparseSpace , class TDenseSpace >
Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::ResidualBasedBossakDisplacementScheme ( const double  Alpha = 0.0)
inlineexplicit

Constructor from a Bossak parameter.

Parameters
Alphais the Bossak parameter. Default value is 0, which is the Newmark method.
Note
The Newmark beta parameter is set to 0.25, for mean constant acceleration.

◆ ResidualBasedBossakDisplacementScheme() [3/4]

template<class TSparseSpace , class TDenseSpace >
Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::ResidualBasedBossakDisplacementScheme ( const double  Alpha,
const double  NewmarkBeta 
)
inlineexplicit

Constructor.

Parameters
AlphaThe Bossak parameter.
NewarkBetaThe Newmark parameter.
Note
A Bossak parameter (Alpha) of 0 results in the Newmark method.
A Newmark parameter (NewmarkBeta) of 0.25 results in mean constant acceleration.

◆ ResidualBasedBossakDisplacementScheme() [4/4]

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

◆ ~ResidualBasedBossakDisplacementScheme()

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

Member Function Documentation

◆ AddDynamicsToLHS()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::AddDynamicsToLHS ( LocalSystemMatrixType LHS_Contribution,
LocalSystemMatrixType D,
LocalSystemMatrixType M,
const ProcessInfo rCurrentProcessInfo 
)
inlineoverrideprotectedvirtual

Add dynamic left hand side contribution from Element s.

\[ M \cdot c_0 + D \cdot c_1 + K \]

Parameters
LHS_ContributionDynamic contribution for the left hand side.
DDamping matrix.
MMass matrix.
rCurrentProcessInfoCurrent ProcessInfo.

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

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

◆ AddDynamicsToRHS() [1/2]

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::AddDynamicsToRHS ( Condition rCondition,
LocalSystemVectorType RHS_Contribution,
LocalSystemMatrixType D,
LocalSystemMatrixType M,
const ProcessInfo rCurrentProcessInfo 
)
inlineoverrideprotectedvirtual

Add dynamic right hand side contribution of a Condition.

\[ b - (1-alpha)*M*a_n+1 - alpha*M*a_n - D*v_n \]

Parameters
rConditionCondition to compute the contribution from.
RHS_ContributionDynamic contribution for the right hand side.
DDamping matrix.
MMass matrix.
rCurrentProcessInfoCurrent ProcessInfo.

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

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

◆ AddDynamicsToRHS() [2/2]

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::AddDynamicsToRHS ( Element rElement,
LocalSystemVectorType RHS_Contribution,
LocalSystemMatrixType D,
LocalSystemMatrixType M,
const ProcessInfo rCurrentProcessInfo 
)
inlineoverrideprotectedvirtual

Add dynamic right hand side contribution from an Element.

\[ b - (1 - \alpha) \cdot M \cdot a_{n+1} - \alpha \cdot M \cdot a_n - D \cdot v_n \]

Parameters
rElementElement to compute the contribution from.
RHS_ContributionDynamic contribution for the right hand side.
DDamping matrix.
MMass matrix.
rCurrentProcessInfoCurrent ProcessInfo.

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

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

◆ AssignSettings()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::AssignSettings ( const Parameters  ThisParameters)
inlineoverrideprotectedvirtual

◆ CalculateBossakCoefficients()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::CalculateBossakCoefficients ( )
inline

Recalculate the Newmark coefficients, taking the alpha parameters into account.

◆ Check()

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

Check whether the scheme and the provided ModelPart are configured correctly.

Checks may be expensive as the function is designed to catch user errors.

Parameters
rModelPartModelPart to check.
Returns
0 if ok, nonzero otherwise.

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

◆ Clear()

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

Release dynamic memory allocated by this instance.

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

◆ Clone()

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

◆ Create()

template<class TSparseSpace , class TDenseSpace >
BaseType::Pointer Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::Create ( Parameters  ThisParameters) const
inlineoverridevirtual

Construct a dynamically allocated new scheme from a Parameters object.

Parameters
ThisParametersThe configuration parameters.

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

Reimplemented in Kratos::ResidualBasedPseudoStaticDisplacementScheme< TSparseSpace, TDenseSpace >, and Kratos::ResidualBasedNewmarkDisplacementScheme< TSparseSpace, TDenseSpace >.

◆ GetDefaultParameters()

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

This function returns the default Parameters to help avoiding conflicts between different constructors.

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

Reimplemented in Kratos::ResidualBasedPseudoStaticDisplacementScheme< TSparseSpace, TDenseSpace >, and Kratos::ResidualBasedNewmarkDisplacementScheme< TSparseSpace, TDenseSpace >.

◆ Info()

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

◆ InitializeSolutionStep()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::InitializeSolutionStep ( ModelPart rModelPart,
TSystemMatrixType A,
TSystemVectorType Dx,
TSystemVectorType b 
)
inlineoverridevirtual

Prepare state variables for a new solution step.

Note
This function should be called before solving a solution step, or restarting a failed one.
Parameters
rModelPartModelPart to update.
ALeft hand side matrix.
DxPrimary variable updates.
bRight hand side vector.

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

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

◆ KRATOS_CLASS_POINTER_DEFINITION()

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

◆ Name()

template<class TSparseSpace , class TDenseSpace >
static std::string Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::Name ( )
inlinestatic

Return the name of the class as used in the settings (snake_case).

◆ Predict()

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

Apply the predictor.

\[ x_{k+1} = x_{k} + v_{k} \cdot \Delta t \]

Parameters
rModelPartModelPart to update.
rDofSetSet of all primary variables.
rALeft hand side matrix.
rDxPrimary variable updates.
rbRight hand side vector.

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

Reimplemented in Kratos::ResidualBasedPseudoStaticDisplacementScheme< TSparseSpace, TDenseSpace >, and Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >.

◆ PrintData()

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

◆ PrintInfo()

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

◆ Update()

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

Update state variables within a newton iteration at the end of the time step.

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

Parameters
rModelPartModelPart to update.
rDofSetSet of all primary variables.
rALeft hand side matrix.
rDxPrimary variable updates.
rbRight hand side vector.

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

Reimplemented in Kratos::ResidualBasedPseudoStaticDisplacementScheme< TSparseSpace, TDenseSpace >, and Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >.

◆ UpdateAcceleration()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::UpdateAcceleration ( array_1d< double, 3 > &  rCurrentAcceleration,
const array_1d< double, 3 > &  rDeltaDisplacement,
const array_1d< double, 3 > &  rPreviousVelocity,
const array_1d< double, 3 > &  rPreviousAcceleration 
)
inlineprotected

Update the second time derivative.

Parameters
rCurrentAccelerationCurrent velocity.
rDeltaDisplacementDisplacement increment.
rPreviousVelocityPrevious velocity.
rPreviousAccelerationPrevious acceleration.

◆ UpdateVelocity()

template<class TSparseSpace , class TDenseSpace >
void Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::UpdateVelocity ( array_1d< double, 3 > &  rCurrentVelocity,
const array_1d< double, 3 > &  rDeltaDisplacement,
const array_1d< double, 3 > &  rPreviousVelocity,
const array_1d< double, 3 > &  rPreviousAcceleration 
)
inlineprotected

Update the first time derivative.

Parameters
rCurrentVelocityCurrent velocity.
rDeltaDisplacementDisplacement increment.
rPreviousVelocityPrevious velocity.
rPreviousAccelerationPrevious acceleration.

Member Data Documentation

◆ mBossak

template<class TSparseSpace , class TDenseSpace >
BossakAlphaMethod Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::mBossak
protected

Bossak Alpha parameters.

◆ mNewmark

template<class TSparseSpace , class TDenseSpace >
NewmarkMethod Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::mNewmark
protected

Newmark Beta parameters.

◆ mpDofUpdater

template<class TSparseSpace , class TDenseSpace >
TSparseSpace::DofUpdaterPointerType Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >::mpDofUpdater = TSparseSpace::CreateDofUpdater()
protected

◆ mVector

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

Aggregate struct for velocities and accelerations.


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