14 #ifndef KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_SLIP_H
15 #define KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_SLIP_H
66 template<
class TSparseSpace,
118 "name" : "ResidualBasedIncrementalUpdateStaticSchemeSlip",
134 unsigned int BlockSize):
145 mpRotationTool(pRotationTool)
166 return Kratos::make_shared<ClassType>(ThisParameters);
178 mpRotationTool->RotateVelocities(r_model_part);
182 mpRotationTool->RecoverVelocities(r_model_part);
199 mpRotationTool->Rotate(LHS_Contribution,RHS_Contribution,rCurrentElement.
GetGeometry());
200 mpRotationTool->ApplySlipCondition(LHS_Contribution,RHS_Contribution,rCurrentElement.
GetGeometry());
215 mpRotationTool->Rotate(RHS_Contribution,rCurrentElement.
GetGeometry());
216 mpRotationTool->ApplySlipCondition(RHS_Contribution,rCurrentElement.
GetGeometry());
232 mpRotationTool->Rotate(LHS_Contribution,
Temp,rCurrentElement.
GetGeometry());
233 mpRotationTool->ApplySlipCondition(LHS_Contribution,
Temp,rCurrentElement.
GetGeometry());
250 mpRotationTool->Rotate(LHS_Contribution,RHS_Contribution,rCurrentCondition.
GetGeometry());
251 mpRotationTool->ApplySlipCondition(LHS_Contribution,RHS_Contribution,rCurrentCondition.
GetGeometry());
266 mpRotationTool->Rotate(RHS_Contribution,rCurrentCondition.
GetGeometry());
267 mpRotationTool->ApplySlipCondition(RHS_Contribution,rCurrentCondition.
GetGeometry());
278 return "static_slip_scheme";
291 std::string
Info()
const override
293 return "ResidualBasedIncrementalUpdateStaticSchemeSlip";
Base class for all Conditions.
Definition: condition.h:59
Base class for all Elements.
Definition: element.h:60
std::vector< std::size_t > EquationIdVectorType
Definition: element.h:98
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
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
int GetInt() const
This method returns the integer contained in the current Parameter.
Definition: kratos_parameters.cpp:666
void ValidateAndAssignDefaults(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing match the form prescribed by th...
Definition: kratos_parameters.cpp:1306
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
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
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Local system vector type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:84
BaseType::TSystemMatrixType TSystemMatrixType
Matrix type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:78
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::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
BaseType::LocalSystemVectorType LocalSystemVectorType
Local system matrix type definition.
Definition: residualbased_incrementalupdate_static_scheme.h:82
Scheme for the solution of problems involving a slip condition.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:70
std::string Info() const override
Turn back information as a string.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:291
CoordinateTransformationUtils< LocalSystemMatrixType, LocalSystemVectorType, double >::Pointer RotationToolPointerType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:95
ResidualBasedIncrementalUpdateStaticSchemeSlip< TSparseSpace, TDenseSpace > ClassType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:82
BaseType::DofsArrayType DofsArrayType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:86
void CalculateRHSContribution(Element &rCurrentElement, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
Obtain an element's local contribution to the RHS and apply slip conditions if needed.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:206
Scheme< TSparseSpace, TDenseSpace > BaseSchemeType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:78
BaseSchemeType::Pointer Create(Parameters ThisParameters) const override
Create method.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:164
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:297
ResidualBasedIncrementalUpdateStaticScheme< TSparseSpace, TDenseSpace > BaseType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:80
void CalculateLHSContribution(Element &rCurrentElement, LocalSystemMatrixType &LHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
Obtain an element's local contribution to the system matrix and apply slip conditions if needed.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:222
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:303
void CalculateRHSContribution(Condition &rCurrentCondition, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
Obtain a condition's local contribution to the RHS and apply slip conditions if needed.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:257
static std::string Name()
Returns the name of the class as used in the settings (snake_case format)
Definition: residualbased_incrementalupdate_static_scheme_slip.h:276
ResidualBasedIncrementalUpdateStaticSchemeSlip(unsigned int DomainSize, unsigned int BlockSize)
Constructor.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:133
void Update(ModelPart &r_model_part, DofsArrayType &rDofSet, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b) override
Update the degrees of freedom after a solution iteration.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:170
BaseType::TSystemVectorType TSystemVectorType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:90
KRATOS_CLASS_POINTER_DEFINITION(ResidualBasedIncrementalUpdateStaticSchemeSlip)
BaseType::TSystemMatrixType TSystemMatrixType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:88
ResidualBasedIncrementalUpdateStaticSchemeSlip(Parameters ThisParameters)
Constructor. The pseudo static scheme (parameters)
Definition: residualbased_incrementalupdate_static_scheme_slip.h:112
BaseType::TDataType TDataType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:84
void CalculateSystemContributions(Condition &rCurrentCondition, LocalSystemMatrixType &LHS_Contribution, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
Obtain a condition's local contribution to the system and apply slip conditions if needed.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:240
BaseType::LocalSystemVectorType LocalSystemVectorType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:92
void CalculateSystemContributions(Element &rCurrentElement, LocalSystemMatrixType &LHS_Contribution, LocalSystemVectorType &RHS_Contribution, Element::EquationIdVectorType &EquationId, const ProcessInfo &CurrentProcessInfo) override
Obtain an element's local contribution to the system and apply slip conditions if needed.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:189
CoordinateTransformationUtils< LocalSystemMatrixType, LocalSystemVectorType, double > RotationToolType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:94
ResidualBasedIncrementalUpdateStaticSchemeSlip()
Default constructor.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:104
~ResidualBasedIncrementalUpdateStaticSchemeSlip() override
Destructor.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:149
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Definition: residualbased_incrementalupdate_static_scheme_slip.h:93
ResidualBasedIncrementalUpdateStaticSchemeSlip(RotationToolPointerType pRotationTool)
Constructor providing a custom rotation tool.
Definition: residualbased_incrementalupdate_static_scheme_slip.h:143
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 void EquationId(const Element &rElement, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo)
This method gets the eqaution id corresponding to the current element.
Definition: scheme.h:636
typename TSparseSpace::VectorType TSystemVectorType
Vector type definition.
Definition: scheme.h:74
typename TDenseSpace::VectorType LocalSystemVectorType
Local system vector type definition.
Definition: scheme.h:80
typename TSparseSpace::DataType TDataType
Data type definition.
Definition: scheme.h:68
typename TDenseSpace::MatrixType LocalSystemMatrixType
Local system matrix type definition.
Definition: scheme.h:77
#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
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
shared_ptr< C > make_shared(Args &&...args)
Definition: smart_pointers.h:40
Temp
Definition: PecletTest.py:105
int domain_size
Definition: face_heat.py:4
b
Definition: generate_total_lagrangian_mixed_volumetric_strain_element.py:31
int block_size
Definition: generate_total_lagrangian_mixed_volumetric_strain_element.py:16
A
Definition: sensitivityMatrix.py:70