13 #if !defined(KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_H )
14 #define KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_H
52 template<
class TSparseSpace,
56 :
public Scheme<TSparseSpace,TDenseSpace>
141 return Kratos::make_shared<ClassType>(ThisParameters);
162 mpDofUpdater->UpdateDofs(rDofSet, rDx);
231 rCurrentCondition.
CalculateLocalSystem(rLHSContribution, rRHSContribution, rCurrentProcessInfo);
312 this->mpDofUpdater->Clear();
323 "name" : "static_scheme"
329 return default_parameters;
338 return "static_scheme";
354 std::string
Info()
const override
356 return "ResidualBasedIncrementalUpdateStaticScheme";
415 typename TSparseSpace::DofUpdaterPointerType mpDofUpdater = TSparseSpace::CreateDofUpdater();
Base class for all Conditions.
Definition: condition.h:59
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: condition.h:260
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:440
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: condition.h:408
Base class for all Elements.
Definition: element.h:60
virtual void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:423
virtual void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:437
virtual void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const
Definition: element.h:258
virtual void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo)
Definition: element.h:405
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
void RecursivelyAddMissingParameters(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing contain at least all parameters...
Definition: kratos_parameters.cpp:1457
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
This class provides the implementation of a static scheme.
Definition: residualbased_incrementalupdate_static_scheme.h:57
void CalculateSystemContributions(Element &rCurrentElement, LocalSystemMatrixType &rLHSContribution, LocalSystemVectorType &rRHSContribution, EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
This function is designed to be called in the builder and solver to introduce the selected time integ...
Definition: residualbased_incrementalupdate_static_scheme.h:196
static std::string Name()
Returns the name of the class as used in the settings (snake_case format)
Definition: residualbased_incrementalupdate_static_scheme.h:336
void CalculateRHSContribution(Element &rCurrentElement, LocalSystemVectorType &rRHSContribution, EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
This function is designed to calculate just the RHS contribution.
Definition: residualbased_incrementalupdate_static_scheme.h:245
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: residualbased_incrementalupdate_static_scheme.h:360
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Local system vector type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:84
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: residualbased_incrementalupdate_static_scheme.h:366
Parameters GetDefaultParameters() const override
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: residualbased_incrementalupdate_static_scheme.h:319
~ResidualBasedIncrementalUpdateStaticScheme() override
Definition: residualbased_incrementalupdate_static_scheme.h:125
BaseType::Pointer Create(Parameters ThisParameters) const override
Create method.
Definition: residualbased_incrementalupdate_static_scheme.h:139
ModelPart::ElementsContainerType ElementsArrayType
Elements containers definition.
Definition: residualbased_incrementalupdate_static_scheme.h:87
ModelPart::ConditionsContainerType ConditionsArrayType
Conditions containers definition.
Definition: residualbased_incrementalupdate_static_scheme.h:89
BaseType::TDataType TDataType
Data type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:76
std::string Info() const override
Turn back information as a string.
Definition: residualbased_incrementalupdate_static_scheme.h:354
BaseType::TSystemMatrixType TSystemMatrixType
Matrix type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:78
Scheme< TSparseSpace, TDenseSpace > BaseType
Base class definition.
Definition: residualbased_incrementalupdate_static_scheme.h:67
KRATOS_CLASS_POINTER_DEFINITION(ResidualBasedIncrementalUpdateStaticScheme)
Pointer definition of ResidualBasedIncrementalUpdateStaticScheme.
Element::EquationIdVectorType EquationIdVectorType
The definition of the vector containing the equation ids.
Definition: residualbased_incrementalupdate_static_scheme.h:92
void Predict(ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
Performing the prediction of the solution.
Definition: residualbased_incrementalupdate_static_scheme.h:174
ResidualBasedIncrementalUpdateStaticScheme(ResidualBasedIncrementalUpdateStaticScheme &rOther)
Definition: residualbased_incrementalupdate_static_scheme.h:118
void Update(ModelPart &rModelPart, DofsArrayType &rDofSet, TSystemMatrixType &rA, TSystemVectorType &rDx, TSystemVectorType &rb) override
Performing the update of the solution.
Definition: residualbased_incrementalupdate_static_scheme.h:152
BaseType::DofsArrayType DofsArrayType
DoF array type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:73
void Clear() override
Liberate internal storage.
Definition: residualbased_incrementalupdate_static_scheme.h:310
void CalculateSystemContributions(Condition &rCurrentCondition, LocalSystemMatrixType &rLHSContribution, LocalSystemVectorType &rRHSContribution, EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: residualbased_incrementalupdate_static_scheme.h:221
void CalculateRHSContribution(Condition &rCurrentCondition, LocalSystemVectorType &rRHSContribution, EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
Functions totally analogous to the precedent but applied to the "condition" objects.
Definition: residualbased_incrementalupdate_static_scheme.h:268
BaseType::TSystemVectorType TSystemVectorType
Vector type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:80
void CalculateLHSContribution(Element &rCurrentElement, LocalSystemMatrixType &rLHSContribution, EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
This function is designed to calculate just the LHS contribution.
Definition: residualbased_incrementalupdate_static_scheme.h:291
ResidualBasedIncrementalUpdateStaticScheme()
Definition: residualbased_incrementalupdate_static_scheme.h:112
ResidualBasedIncrementalUpdateStaticScheme(Parameters ThisParameters)
Constructor. The pseudo static scheme (parameters)
Definition: residualbased_incrementalupdate_static_scheme.h:102
BaseType::LocalSystemVectorType LocalSystemVectorType
Local system matrix type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:82
ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace > ClassType
Definition: residualbased_incrementalupdate_static_scheme.h:70
This class provides the implementation of the basic tasks that are needed by the solution strategy.
Definition: scheme.h:56
typename TSparseSpace::MatrixType TSystemMatrixType
Matrix type definition.
Definition: scheme.h:71
virtual Parameters GetDefaultParameters() const
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: scheme.h:693
typename TSparseSpace::VectorType TSystemVectorType
Vector type definition.
Definition: scheme.h:74
typename TDenseSpace::VectorType LocalSystemVectorType
Local system vector type definition.
Definition: scheme.h:80
virtual Parameters ValidateAndAssignParameters(Parameters ThisParameters, const Parameters DefaultParameters) const
This method validate and assign default parameters.
Definition: scheme.h:773
typename TSparseSpace::DataType TDataType
Data type definition.
Definition: scheme.h:68
typename TDenseSpace::MatrixType LocalSystemMatrixType
Local system matrix type definition.
Definition: scheme.h:77
virtual void AssignSettings(const Parameters ThisParameters)
This method assigns settings to member variables.
Definition: scheme.h:786
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21