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

This class provides the implementation of the basic tasks that are needed by the solution strategy. More...

#include <scheme.h>

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

Public Member Functions

Life Cycle
 Scheme ()
 Default Constructor. More...
 
 Scheme (Parameters ThisParameters)
 Constructor with Parameters. More...
 
 Scheme (Scheme &rOther)
 
virtual ~Scheme ()
 
Input and output
virtual std::string Info () const
 Turn back information as a string. More...
 
virtual void PrintInfo (std::ostream &rOStream) const
 Print information about this object. More...
 
virtual void PrintData (std::ostream &rOStream) const
 Print object's data. More...
 

Protected Member Functions

Protected Operations
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
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 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...
 
 KRATOS_CLASS_POINTER_DEFINITION (Scheme)
 Pointer definition of Scheme. More...
 

Operations

virtual ClassType::Pointer Create (Parameters ThisParameters) const
 Create method. More...
 
virtual Pointer Clone ()
 Clone 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 InitializeSolutionStep (ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
 Function called once at the beginning of each solution step. 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 Predict (ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
 Performing the prediction of the solution. More...
 
virtual void Update (ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
 Performing the update of the solution. 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 void Clear ()
 Liberate internal storage. More...
 
virtual int Check (const ModelPart &rModelPart) const
 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. More...
 
virtual int Check (ModelPart &rModelPart)
 
virtual void CalculateSystemContributions (Element &rElement, LocalSystemMatrixType &LHS_Contribution, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
 This function is designed to be called in the builder and solver to introduce the selected time integration scheme. More...
 
virtual void CalculateSystemContributions (Condition &rCondition, LocalSystemMatrixType &LHS_Contribution, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
 Functions totally analogous to the precedent but applied to the "condition" objects. More...
 
virtual void CalculateRHSContribution (Element &rElement, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
 This function is designed to calculate just the RHS contribution. More...
 
virtual void CalculateRHSContribution (Condition &rCondition, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &rEquationIdVector, const ProcessInfo &rCurrentProcessInfo)
 Functions totally analogous to the precedent but applied to the "condition" objects. More...
 
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...
 
virtual Parameters GetDefaultParameters () const
 This method provides the defaults parameters to avoid conflicts between the different constructors. More...
 
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::Scheme< TSparseSpace, TDenseSpace >

This class provides the implementation of the basic tasks that are needed by the solution strategy.

It is intended to be the place for tailoring the solution strategies to problem specific tasks.

Template Parameters
TSparseSpaceThe sparse space considered
TDenseSpaceThe dense space considered
Author
Riccardo Rossi

Member Typedef Documentation

◆ ClassType

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

The definition of the current class.

◆ ConditionsArrayType

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

Conditions containers definition.

◆ DofsArrayType

template<class TSparseSpace , class TDenseSpace >
using Kratos::Scheme< TSparseSpace, TDenseSpace >::DofsArrayType = ModelPart::DofsArrayType

DoF array type definition.

◆ ElementsArrayType

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

Elements containers definition.

◆ LocalSystemMatrixType

template<class TSparseSpace , class TDenseSpace >
using Kratos::Scheme< TSparseSpace, TDenseSpace >::LocalSystemMatrixType = typename TDenseSpace::MatrixType

Local system matrix type definition.

◆ LocalSystemVectorType

template<class TSparseSpace , class TDenseSpace >
using Kratos::Scheme< TSparseSpace, TDenseSpace >::LocalSystemVectorType = typename TDenseSpace::VectorType

Local system vector type definition.

◆ TDataType

template<class TSparseSpace , class TDenseSpace >
using Kratos::Scheme< TSparseSpace, TDenseSpace >::TDataType = typename TSparseSpace::DataType

Data type definition.

◆ TDofType

template<class TSparseSpace , class TDenseSpace >
using Kratos::Scheme< TSparseSpace, TDenseSpace >::TDofType = Dof<double>

DoF type definition.

◆ TSystemMatrixType

template<class TSparseSpace , class TDenseSpace >
using Kratos::Scheme< TSparseSpace, TDenseSpace >::TSystemMatrixType = typename TSparseSpace::MatrixType

Matrix type definition.

◆ TSystemVectorType

template<class TSparseSpace , class TDenseSpace >
using Kratos::Scheme< TSparseSpace, TDenseSpace >::TSystemVectorType = typename TSparseSpace::VectorType

Vector type definition.

Constructor & Destructor Documentation

◆ Scheme() [1/3]

template<class TSparseSpace , class TDenseSpace >
Kratos::Scheme< TSparseSpace, TDenseSpace >::Scheme ( )
inlineexplicit

Default Constructor.

Initializes the flags

◆ Scheme() [2/3]

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

Constructor with Parameters.

◆ Scheme() [3/3]

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

Copy Constructor.

◆ ~Scheme()

template<class TSparseSpace , class TDenseSpace >
virtual Kratos::Scheme< TSparseSpace, TDenseSpace >::~Scheme ( )
inlinevirtual

Destructor.

Member Function Documentation

◆ AssignSettings()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::AssignSettings ( const Parameters  ThisParameters)
inlineprotectedvirtual

◆ CalculateLHSContribution() [1/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::CalculateLHSContribution ( Condition rCondition,
LocalSystemMatrixType LHS_Contribution,
Element::EquationIdVectorType rEquationIdVector,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

◆ CalculateLHSContribution() [2/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::CalculateLHSContribution ( Element rElement,
LocalSystemMatrixType LHS_Contribution,
Element::EquationIdVectorType rEquationIdVector,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

◆ CalculateOutputData()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::CalculateOutputData ( ModelPart rModelPart,
DofsArrayType rDofSet,
TSystemMatrixType A,
TSystemVectorType Dx,
TSystemVectorType b 
)
inlinevirtual

Functions to be called to prepare the data needed for the output of results.

Warning
Must be defined in derived classes
Parameters
rModelPartThe model part of the problem to solve
rDofSetSet of all primary variables
ALHS matrix
DxIncremental update of primary variables
bRHS Vector

◆ CalculateRHSContribution() [1/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::CalculateRHSContribution ( Condition rCondition,
LocalSystemVectorType RHS_Contribution,
Element::EquationIdVectorType rEquationIdVector,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

Functions totally analogous to the precedent but applied to the "condition" objects.

Parameters
rConditionThe condition to compute
RHS_ContributionThe RHS vector contribution
rEquationIdVectorThe ID's of the condition degrees of freedom
rCurrentProcessInfoThe current process info instance

Reimplemented in Kratos::ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdamsMoultonScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointStaticScheme< TSparseSpace, TDenseSpace >, Kratos::NewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::FluxCorrectedShallowWaterScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedSimpleSteadyScheme< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitMultiStageKimScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitCentralDifferencesScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitCDScheme< TSparseSpace, TDenseSpace >, Kratos::PFEM2MonolithicSlipScheme< TSparseSpace, TDenseSpace >, Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::MPMExplicitScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorVelocityBossakSchemeTurbulent< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedRelaxationScheme< TSparseSpace, TDenseSpace >, Kratos::EigensolverDynamicScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, Kratos::GeneralizedNewmarkGN11Scheme< TSparseSpace, TDenseSpace >, Kratos::SteadyScalarScheme< TSparseSpace, TDenseSpace >, and Kratos::AlgebraicFluxCorrectedSteadyScalarScheme< TSparseSpace, TDenseSpace >.

◆ CalculateRHSContribution() [2/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::CalculateRHSContribution ( Element rElement,
LocalSystemVectorType RHS_Contribution,
Element::EquationIdVectorType rEquationIdVector,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

This function is designed to calculate just the RHS contribution.

Parameters
rElementThe element to compute
RHS_ContributionThe RHS vector contribution
rEquationIdVectorThe ID's of the element degrees of freedom
rCurrentProcessInfoThe current process info instance

Reimplemented in Kratos::SteadyScalarScheme< TSparseSpace, TDenseSpace >, Kratos::AlgebraicFluxCorrectedSteadyScalarScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdamsMoultonScheme< TSparseSpace, TDenseSpace >, Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::FluxCorrectedShallowWaterScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedSimpleSteadyScheme< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitMultiStageKimScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitCentralDifferencesScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitVVScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitCDScheme< TSparseSpace, TDenseSpace >, Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::MPMExplicitScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedRelaxationScheme< TSparseSpace, TDenseSpace >, Kratos::EigensolverDynamicScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticDampedUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PFEM2MonolithicSlipScheme< TSparseSpace, TDenseSpace >, Kratos::NewmarkQuasistaticDampedUPwScheme< TSparseSpace, TDenseSpace >, Kratos::NewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, Kratos::GeneralizedNewmarkGN11Scheme< TSparseSpace, TDenseSpace >, and Kratos::ResidualBasedPredictorCorrectorVelocityBossakSchemeTurbulent< TSparseSpace, TDenseSpace >.

◆ CalculateSystemContributions() [1/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::CalculateSystemContributions ( Condition rCondition,
LocalSystemMatrixType LHS_Contribution,
LocalSystemVectorType RHS_Contribution,
Element::EquationIdVectorType rEquationIdVector,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

Functions totally analogous to the precedent but applied to the "condition" objects.

Parameters
rConditionThe condition to compute
LHS_ContributionThe LHS matrix contribution
RHS_ContributionThe RHS vector contribution
rEquationIdVectorThe ID's of the condition degrees of freedom
rCurrentProcessInfoThe current process info instance

Reimplemented in Kratos::ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdamsMoultonScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::FluxCorrectedShallowWaterScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointSteadyScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointBossakScheme< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentScheme< TSparseSpace, TDenseSpace >, Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedRelaxationScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PFEM2MonolithicSlipScheme< TSparseSpace, TDenseSpace >, Kratos::NewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, Kratos::GeneralizedNewmarkGN11Scheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorVelocityBossakSchemeTurbulent< TSparseSpace, TDenseSpace >, Kratos::DamUPScheme< TSparseSpace, TDenseSpace >, Kratos::DamPScheme< TSparseSpace, TDenseSpace >, Kratos::EigensolverDynamicScheme< TSparseSpace, TDenseSpace >, Kratos::EigensolverNitscheStabilizationScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedSimpleSteadyScheme< TSparseSpace, TDenseSpace >, Kratos::AlgebraicFluxCorrectedSteadyScalarScheme< TSparseSpace, TDenseSpace >, and Kratos::SteadyScalarScheme< TSparseSpace, TDenseSpace >.

◆ CalculateSystemContributions() [2/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::CalculateSystemContributions ( Element rElement,
LocalSystemMatrixType LHS_Contribution,
LocalSystemVectorType RHS_Contribution,
Element::EquationIdVectorType rEquationIdVector,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

This function is designed to be called in the builder and solver to introduce the selected time integration scheme.

It "asks" the matrix needed to the element and performs the operations needed to introduce the selected time integration scheme. This function calculates at the same time the contribution to the LHS and to the RHS of the system

Parameters
rElementThe element to compute
LHS_ContributionThe LHS matrix contribution
RHS_ContributionThe RHS vector contribution
rEquationIdVectorThe ID's of the element degrees of freedom
rCurrentProcessInfoThe current process info instance

Reimplemented in Kratos::SteadyScalarScheme< TSparseSpace, TDenseSpace >, Kratos::AlgebraicFluxCorrectedSteadyScalarScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedIncrementalUpdateStaticSIMPScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdamsMoultonScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointSteadyScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::FluxCorrectedShallowWaterScheme< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentScheme< TSparseSpace, TDenseSpace >, Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedRelaxationScheme< TSparseSpace, TDenseSpace >, Kratos::EigensolverDynamicScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticDampedUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PFEM2MonolithicSlipScheme< TSparseSpace, TDenseSpace >, Kratos::EigensolverNitscheStabilizationScheme< TSparseSpace, TDenseSpace >, Kratos::NewmarkQuasistaticDampedUPwScheme< TSparseSpace, TDenseSpace >, Kratos::NewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, Kratos::GeneralizedNewmarkGN11Scheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedSimpleSteadyScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorVelocityBossakSchemeTurbulent< TSparseSpace, TDenseSpace >, Kratos::IncrementalUpdateStaticDampedSmoothingScheme< TSparseSpace, TDenseSpace >, Kratos::DamUPScheme< TSparseSpace, TDenseSpace >, and Kratos::DamPScheme< TSparseSpace, TDenseSpace >.

◆ Check() [1/2]

template<class TSparseSpace , class TDenseSpace >
virtual int Kratos::Scheme< TSparseSpace, TDenseSpace >::Check ( const ModelPart rModelPart) const
inlinevirtual

◆ Check() [2/2]

template<class TSparseSpace , class TDenseSpace >
virtual int Kratos::Scheme< TSparseSpace, TDenseSpace >::Check ( ModelPart rModelPart)
inlinevirtual

◆ Clean()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::Clean ( )
inlinevirtual

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.

Warning
Must be implemented in the derived classes

◆ CleanOutputData()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::CleanOutputData ( )
inlinevirtual

Functions that cleans the results data.

Warning
Must be implemented in the derived classes

◆ Clear()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::Clear ( )
inlinevirtual

◆ Clone()

template<class TSparseSpace , class TDenseSpace >
virtual Pointer Kratos::Scheme< TSparseSpace, TDenseSpace >::Clone ( )
inlinevirtual

◆ ConditionsAreInitialized()

template<class TSparseSpace , class TDenseSpace >
bool Kratos::Scheme< TSparseSpace, TDenseSpace >::ConditionsAreInitialized ( )
inline

This method returns if the conditions are initialized.

Returns
True if initialized, false otherwise

◆ Create()

template<class TSparseSpace , class TDenseSpace >
virtual ClassType::Pointer Kratos::Scheme< TSparseSpace, TDenseSpace >::Create ( Parameters  ThisParameters) const
inlinevirtual

◆ ElementsAreInitialized()

template<class TSparseSpace , class TDenseSpace >
bool Kratos::Scheme< TSparseSpace, TDenseSpace >::ElementsAreInitialized ( )
inline

This method returns if the elements are initialized.

Returns
True if initialized, false otherwise

◆ EquationId() [1/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::EquationId ( const Condition rCondition,
Element::EquationIdVectorType rEquationId,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

Functions totally analogous to the precedent but applied to the "condition" objects.

Parameters
rConditionThe condition to compute
rEquationIdThe ID's of the condition degrees of freedom
rCurrentProcessInfoThe current process info instance

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

◆ EquationId() [2/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::EquationId ( const Element rElement,
Element::EquationIdVectorType rEquationId,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

This method gets the eqaution id corresponding to the current element.

Parameters
rElementThe element to compute
rEquationIdThe ID's of the element degrees of freedom
rCurrentProcessInfoThe current process info instance

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

◆ FinalizeNonLinIteration()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::FinalizeNonLinIteration ( ModelPart rModelPart,
TSystemMatrixType A,
TSystemVectorType Dx,
TSystemVectorType b 
)
inlinevirtual

◆ FinalizeSolutionStep()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::FinalizeSolutionStep ( ModelPart rModelPart,
TSystemMatrixType A,
TSystemVectorType Dx,
TSystemVectorType b 
)
inlinevirtual

◆ GetDefaultParameters()

template<class TSparseSpace , class TDenseSpace >
virtual Parameters Kratos::Scheme< TSparseSpace, TDenseSpace >::GetDefaultParameters ( ) const
inlinevirtual

◆ GetDofList() [1/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::GetDofList ( const Condition rCondition,
Element::DofsVectorType rDofList,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

Function that returns the list of Degrees of freedom to be assembled in the system for a Given condition.

Parameters
rConditionThe condition to compute
rDofListThe list containing the condition degrees of freedom
rCurrentProcessInfoThe current process info instance

Reimplemented in Kratos::MPMExplicitScheme< TSparseSpace, TDenseSpace >, and Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >.

◆ GetDofList() [2/2]

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::GetDofList ( const Element rElement,
Element::DofsVectorType rDofList,
const ProcessInfo rCurrentProcessInfo 
)
inlinevirtual

Function that returns the list of Degrees of freedom to be assembled in the system for a Given element.

Parameters
pCurrentElementThe element to compute
rDofListThe list containing the element degrees of freedom
rCurrentProcessInfoThe current process info instance

Reimplemented in Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, and Kratos::MPMExplicitScheme< TSparseSpace, TDenseSpace >.

◆ Info()

template<class TSparseSpace , class TDenseSpace >
virtual std::string Kratos::Scheme< TSparseSpace, TDenseSpace >::Info ( ) const
inlinevirtual

◆ Initialize()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::Initialize ( ModelPart rModelPart)
inlinevirtual

This is the place to initialize the Scheme.

This is intended to be called just once when the strategy is initialized

Parameters
rModelPartThe model part of the problem to solve

Reimplemented in Kratos::ResidualBasedBDFCustomScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitMultiStageKimScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitCentralDifferencesScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitCentralDifferencesScheme< TSparseSpace, TDenseSpace >, Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdamsMoultonScheme< TSparseSpace, TDenseSpace >, Kratos::FluxCorrectedShallowWaterScheme< TSparseSpace, TDenseSpace >, Kratos::AlgebraicFluxCorrectedSteadyScalarScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitCDScheme< TSparseSpace, TDenseSpace >, Kratos::MPMExplicitScheme< TSparseSpace, TDenseSpace >, Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, Kratos::VelocityBossakAdjointScheme< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::GeneralizedNewmarkGN11Scheme< TSparseSpace, TDenseSpace >, Kratos::IncrementalUpdateStaticDampedSmoothingScheme< TSparseSpace, TDenseSpace >, Kratos::DamUPScheme< TSparseSpace, TDenseSpace >, Kratos::DamPScheme< TSparseSpace, TDenseSpace >, Kratos::BossakDisplacementSmoothingScheme< TSparseSpace, TDenseSpace >, and Kratos::ExplicitHamiltonScheme< TSparseSpace, TDenseSpace >.

◆ InitializeConditions()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::InitializeConditions ( ModelPart rModelPart)
inlinevirtual

This is the place to initialize the conditions.

This is intended to be called just once when the strategy is initialized

Parameters
rModelPartThe model part of the problem to solve

Reimplemented in Kratos::ResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, and Kratos::ResidualBasedBossakDisplacementRotationScheme< TSparseSpace, TDenseSpace >.

◆ InitializeElements()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::InitializeElements ( ModelPart rModelPart)
inlinevirtual

This is the place to initialize the elements.

This is intended to be called just once when the strategy is initialized

Parameters
rModelPartThe model part of the problem to solve

Reimplemented in Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitCDScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, and Kratos::ResidualBasedBossakDisplacementRotationScheme< TSparseSpace, TDenseSpace >.

◆ InitializeNonLinIteration()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::InitializeNonLinIteration ( ModelPart rModelPart,
TSystemMatrixType A,
TSystemVectorType Dx,
TSystemVectorType b 
)
inlinevirtual

◆ InitializeSolutionStep()

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

Function called once at the beginning of each solution step.

The basic operations to be carried in there are the following:

  • managing variables to be kept constant over the time step (for example time-Scheme constants depending on the actual time step)
    Parameters
    rModelPartThe model part of the problem to solve
    ALHS matrix
    DxIncremental update of primary variables
    bRHS Vector

Reimplemented in Kratos::ResidualBasedImplicitTimeScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBDFDisplacementScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBDFCustomScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitMultiStageKimScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitCentralDifferencesScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdamsMoultonScheme< TSparseSpace, TDenseSpace >, Kratos::FluxCorrectedShallowWaterScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitCDScheme< TSparseSpace, TDenseSpace >, Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentSchemeDEMCoupled< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentScheme< TSparseSpace, TDenseSpace >, Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointBossakScheme< TSparseSpace, TDenseSpace >, Kratos::TrilinosResidualBasedIncrementalAitkenStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorVelocityBossakSchemeTurbulentDEMCoupled< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedRelaxationScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::PFEM2MonolithicSlipScheme< TSparseSpace, TDenseSpace >, Kratos::MPMExplicitScheme< TSparseSpace, TDenseSpace >, Kratos::GeneralizedNewmarkGN11Scheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorVelocityBossakSchemeTurbulent< TSparseSpace, TDenseSpace >, Kratos::DamUPScheme< TSparseSpace, TDenseSpace >, Kratos::DamPScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitHamiltonScheme< TSparseSpace, TDenseSpace >, and Kratos::ResidualBasedBossakDisplacementRotationScheme< TSparseSpace, TDenseSpace >.

◆ KRATOS_CLASS_POINTER_DEFINITION()

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

Pointer definition of Scheme.

◆ Name()

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

Returns the name of the class as used in the settings (snake_case format)

Returns
The name of the class

◆ Predict()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::Predict ( ModelPart rModelPart,
DofsArrayType rDofSet,
TSystemMatrixType A,
TSystemVectorType Dx,
TSystemVectorType b 
)
inlinevirtual

Performing the prediction of the solution.

Warning
Must be defined in derived classes
Parameters
rModelPartThe model part of the problem to solve
ALHS matrix
DxIncremental update of primary variables
bRHS Vector

Reimplemented in Kratos::ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPseudoStaticDisplacementScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdamsMoultonScheme< TSparseSpace, TDenseSpace >, Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBDFDisplacementScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBDFCustomScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitMultiStageKimScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitCentralDifferencesScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitVVScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitCDScheme< TSparseSpace, TDenseSpace >, Kratos::NewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentScheme< TSparseSpace, TDenseSpace >, Kratos::IncrementalUpdateStaticDampedSmoothingScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorVelocityBossakAleScheme< TSparseSpace, TDenseSpace >, Kratos::BackwardEulerMonolithicAleScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorVelocityBossakSchemeTurbulent< TSparseSpace, TDenseSpace >, Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedRelaxationScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkDynamicUPwScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::GeneralizedNewmarkGN11Scheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitForwardEulerScheme< TSparseSpace, TDenseSpace >, Kratos::DamUPScheme< TSparseSpace, TDenseSpace >, Kratos::DamPScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorBossakScheme< TSparseSpace, TDenseSpace >, and Kratos::ResidualBasedBossakDisplacementRotationScheme< TSparseSpace, TDenseSpace >.

◆ PrintData()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::PrintData ( std::ostream &  rOStream) const
inlinevirtual

◆ PrintInfo()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::PrintInfo ( std::ostream &  rOStream) const
inlinevirtual

◆ SchemeIsInitialized()

template<class TSparseSpace , class TDenseSpace >
bool Kratos::Scheme< TSparseSpace, TDenseSpace >::SchemeIsInitialized ( )
inline

This method returns if the scheme is initialized.

Returns
True if initialized, false otherwise

◆ SetConditionsAreInitialized()

template<class TSparseSpace , class TDenseSpace >
void Kratos::Scheme< TSparseSpace, TDenseSpace >::SetConditionsAreInitialized ( bool  ConditionsAreInitializedFlag = true)
inline

This method sets if the conditions have been initialized or not (true by default)

Parameters
ConditionsAreInitializedFlagIf the flag must be set to true or false

◆ SetElementsAreInitialized()

template<class TSparseSpace , class TDenseSpace >
void Kratos::Scheme< TSparseSpace, TDenseSpace >::SetElementsAreInitialized ( bool  ElementsAreInitializedFlag = true)
inline

This method sets if the elements have been initialized or not (true by default)

Parameters
ElementsAreInitializedFlagIf the flag must be set to true or false

◆ SetSchemeIsInitialized()

template<class TSparseSpace , class TDenseSpace >
void Kratos::Scheme< TSparseSpace, TDenseSpace >::SetSchemeIsInitialized ( bool  SchemeIsInitializedFlag = true)
inline

This method sets if the elements have been initialized or not (true by default)

Parameters
ElementsAreInitializedFlagIf the flag must be set to true or false

◆ Update()

template<class TSparseSpace , class TDenseSpace >
virtual void Kratos::Scheme< TSparseSpace, TDenseSpace >::Update ( ModelPart rModelPart,
DofsArrayType rDofSet,
TSystemMatrixType A,
TSystemVectorType Dx,
TSystemVectorType b 
)
inlinevirtual

Performing the update of the solution.

Warning
Must be defined in derived classes
Parameters
rModelPartThe model part of the problem to solve
rDofSetSet of all primary variables
ALHS matrix
DxIncremental update of primary variables
bRHS Vector

Reimplemented in Kratos::ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPseudoStaticDisplacementScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBossakDisplacementScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitMultiStageKimScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitCentralDifferencesScheme< TSparseSpace, TDenseSpace >, Kratos::ShallowWaterResidualBasedBDFScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdamsMoultonScheme< TSparseSpace, TDenseSpace >, Kratos::SteadyScalarScheme< TSparseSpace, TDenseSpace >, Kratos::AlgebraicFluxCorrectedSteadyScalarScheme< TSparseSpace, TDenseSpace >, Kratos::PoroExplicitCDScheme< TSparseSpace, TDenseSpace >, Kratos::MPMResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedSimpleSteadyScheme< TSparseSpace, TDenseSpace >, Kratos::BDF2TurbulentScheme< TSparseSpace, TDenseSpace >, Kratos::IncrementalUpdateStaticDampedSmoothingScheme< TSparseSpace, TDenseSpace >, Kratos::BackwardEulerMonolithicAleScheme< TSparseSpace, TDenseSpace >, Kratos::GeoMechanicsTimeIntegrationScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedAdjointBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace >, Kratos::TrilinosResidualBasedIncrementalAitkenStaticScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedRelaxationScheme< TSparseSpace, TDenseSpace >, Kratos::PoroNewmarkQuasistaticUPwScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBossakScheme< TSparseSpace, TDenseSpace >, Kratos::MPMExplicitScheme< TSparseSpace, TDenseSpace >, Kratos::GeneralizedNewmarkGN11Scheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitForwardEulerScheme< TSparseSpace, TDenseSpace >, Kratos::DamUPScheme< TSparseSpace, TDenseSpace >, Kratos::DamPScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorBossakScheme< TSparseSpace, TDenseSpace >, Kratos::ExplicitHamiltonScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedBossakDisplacementRotationScheme< TSparseSpace, TDenseSpace >, Kratos::ResidualBasedPredictorCorrectorVelocityBossakAleScheme< TSparseSpace, TDenseSpace >, Kratos::PFEM2MonolithicSlipScheme< TSparseSpace, TDenseSpace >, and Kratos::ResidualBasedPredictorCorrectorVelocityBossakSchemeTurbulent< TSparseSpace, TDenseSpace >.

◆ ValidateAndAssignParameters()

template<class TSparseSpace , class TDenseSpace >
virtual Parameters Kratos::Scheme< TSparseSpace, TDenseSpace >::ValidateAndAssignParameters ( Parameters  ThisParameters,
const Parameters  DefaultParameters 
) const
inlineprotectedvirtual

This method validate and assign default parameters.

Parameters
rParametersParameters to be validated
DefaultParametersThe default parameters
Returns
Returns validated Parameters

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

Member Data Documentation

◆ mConditionsAreInitialized

template<class TSparseSpace , class TDenseSpace >
bool Kratos::Scheme< TSparseSpace, TDenseSpace >::mConditionsAreInitialized
protected

Flag taking in account if the elements were initialized correctly or not.

◆ mElementsAreInitialized

template<class TSparseSpace , class TDenseSpace >
bool Kratos::Scheme< TSparseSpace, TDenseSpace >::mElementsAreInitialized
protected

Flag to be used in controlling if the Scheme has been initialized or not.

◆ mSchemeIsInitialized

template<class TSparseSpace , class TDenseSpace >
bool Kratos::Scheme< TSparseSpace, TDenseSpace >::mSchemeIsInitialized
protected

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