14 #if !defined(KRATOS_RESIDUAL_BASED_IMPLICIT_TIME_SCHEME )
15 #define KRATOS_RESIDUAL_BASED_IMPLICIT_TIME_SCHEME
52 template<
class TSparseSpace,
class TDenseSpace >
54 :
public Scheme<TSparseSpace,TDenseSpace>
131 typename BaseType::Pointer
Clone()
override
133 return Kratos::make_shared<ResidualBasedImplicitTimeScheme>(*
this);
168 TCalculateSystemContributions(rCurrentElement, rLHSContribution, rRHSContribution, rEquationId, rCurrentProcessInfo);
170 KRATOS_CATCH(
"ResidualBasedImplicitTimeScheme.CalculateSystemContributions");
189 TCalculateRHSContribution(rCurrentElement, rRHSContribution, rEquationId, rCurrentProcessInfo);
191 KRATOS_CATCH(
"ResidualBasedImplicitTimeScheme.CalculateRHSContribution");
212 TCalculateSystemContributions(rCurrentCondition, rLHSContribution, rRHSContribution, rEquationId, rCurrentProcessInfo);
214 KRATOS_CATCH(
"ResidualBasedImplicitTimeScheme.CalculateSystemContributions");
233 TCalculateRHSContribution(rCurrentCondition, rRHSContribution, rEquationId, rCurrentProcessInfo);
235 KRATOS_CATCH(
"ResidualBasedImplicitTimeScheme.CalculateRHSContribution");
258 const double delta_time = r_current_process_info[DELTA_TIME];
260 KRATOS_ERROR_IF(
delta_time < 1.0e-24) <<
"ERROR:: Detected delta_time = 0 in the Solution Scheme DELTA_TIME. PLEASE : check if the time step is created correctly for the current time step" << std::endl;
262 KRATOS_CATCH(
"ResidualBasedImplicitTimeScheme.InitializeSolutionStep");
278 if(err!=0)
return err;
293 "name" : "residualbased_implicit_time_scheme"
299 return default_parameters;
315 std::string
Info()
const override
317 return "ResidualBasedImplicitTimeScheme";
347 std::vector< Matrix >
M;
348 std::vector< Matrix >
D;
375 KRATOS_ERROR <<
"YOU ARE CALLING THE BASE CLASS OF AddDynamicsToLHS" << std::endl;
395 KRATOS_ERROR <<
"YOU ARE CALLING THE BASE CLASS OF AddDynamicsToRHS" << std::endl;
414 KRATOS_ERROR <<
"YOU ARE CALLING THE BASE CLASS OF AddDynamicsToRHS" << std::endl;
454 template <
class TObjectType>
455 void TCalculateSystemContributions(
456 TObjectType& rObject,
467 rObject.CalculateLocalSystem(rLHSContribution,rRHSContribution,rCurrentProcessInfo);
469 rObject.EquationIdVector(
EquationId,rCurrentProcessInfo);
471 rObject.CalculateMassMatrix(
mMatrix.
M[this_thread],rCurrentProcessInfo);
473 rObject.CalculateDampingMatrix(
mMatrix.
D[this_thread],rCurrentProcessInfo);
479 KRATOS_CATCH(
"ResidualBasedImplicitTimeScheme.TCalculateSystemContributions");
489 template <
class TObjectType>
490 void TCalculateRHSContribution(
491 TObjectType& rObject,
494 const ProcessInfo& rCurrentProcessInfo
501 rObject.CalculateRightHandSide(rRHSContribution,rCurrentProcessInfo);
503 rObject.CalculateMassMatrix(
mMatrix.
M[this_thread], rCurrentProcessInfo);
505 rObject.CalculateDampingMatrix(
mMatrix.
D[this_thread],rCurrentProcessInfo);
507 rObject.EquationIdVector(rEquationId,rCurrentProcessInfo);
511 KRATOS_CATCH(
"ResidualBasedImplicitTimeScheme.TCalculateRHSContribution");
Base class for all Conditions.
Definition: condition.h:59
Base class for all Elements.
Definition: element.h:60
std::vector< DofType::Pointer > DofsVectorType
Definition: element.h:100
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
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
static int ThisThread()
Wrapper for omp_get_thread_num().
Definition: openmp_utils.h:108
static int GetNumThreads()
Returns the current number of threads.
Definition: parallel_utilities.cpp:34
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 is the base class for the implicit time schemes.
Definition: residual_based_implicit_time_scheme.h:55
void CalculateRHSContribution(Condition &rCurrentCondition, LocalSystemVectorType &rRHSContribution, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
Functions that calculates the RHS of a "condition" object.
Definition: residual_based_implicit_time_scheme.h:224
Parameters GetDefaultParameters() const override
This method provides the defaults parameters to avoid conflicts between the different constructors.
Definition: residual_based_implicit_time_scheme.h:289
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.
Definition: residual_based_implicit_time_scheme.h:273
void CalculateRHSContribution(Element &rCurrentElement, LocalSystemVectorType &rRHSContribution, Element::EquationIdVectorType &rEquationId, const ProcessInfo &rCurrentProcessInfo) override
This function is designed to calculate just the RHS contribution.
Definition: residual_based_implicit_time_scheme.h:180
ResidualBasedImplicitTimeScheme()
Definition: residual_based_implicit_time_scheme.h:99
BaseType::Pointer Clone() override
Definition: residual_based_implicit_time_scheme.h:131
ModelPart::ConditionsContainerType ConditionsArrayType
Conditions containers definition.
Definition: residual_based_implicit_time_scheme.h:86
virtual void AddDynamicsToRHS(Element &rCurrentElement, LocalSystemVectorType &rRHSContribution, LocalSystemMatrixType &D, LocalSystemMatrixType &M, const ProcessInfo &rCurrentProcessInfo)
It adds the dynamic RHS contribution of the elements b - M*a - D*v.
Definition: residual_based_implicit_time_scheme.h:387
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.
Definition: residual_based_implicit_time_scheme.h:245
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: residual_based_implicit_time_scheme.h:327
ResidualBasedImplicitTimeScheme(ResidualBasedImplicitTimeScheme &rOther)
Definition: residual_based_implicit_time_scheme.h:122
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 integ...
Definition: residual_based_implicit_time_scheme.h:158
ResidualBasedImplicitTimeScheme(Parameters ThisParameters)
Constructor. The implicit method method.
Definition: residual_based_implicit_time_scheme.h:113
BaseType::TDataType TDataType
Data type definition.
Definition: residual_based_implicit_time_scheme.h:71
ModelPart::NodesContainerType NodesArrayType
Nodes containers definition.
Definition: residual_based_implicit_time_scheme.h:82
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: residual_based_implicit_time_scheme.h:321
GeneralMatrices mMatrix
Definition: residual_based_implicit_time_scheme.h:351
Scheme< TSparseSpace, TDenseSpace > BaseType
Base class definition.
Definition: residual_based_implicit_time_scheme.h:63
virtual void AddDynamicsToLHS(LocalSystemMatrixType &rLHSContribution, LocalSystemMatrixType &D, LocalSystemMatrixType &M, const ProcessInfo &rCurrentProcessInfo)
It adds the dynamic LHS contribution of the elements LHS = d(-RHS)/d(un0) = c0*c0*M + c0*D + K.
Definition: residual_based_implicit_time_scheme.h:368
~ResidualBasedImplicitTimeScheme() override
Definition: residual_based_implicit_time_scheme.h:139
BaseType::TSystemMatrixType TSystemMatrixType
Matrix type definition.
Definition: residual_based_implicit_time_scheme.h:73
BaseType::LocalSystemMatrixType LocalSystemMatrixType
Local system vector type definition.
Definition: residual_based_implicit_time_scheme.h:79
std::string Info() const override
Turn back information as a string.
Definition: residual_based_implicit_time_scheme.h:315
BaseType::LocalSystemVectorType LocalSystemVectorType
Local system matrix type definition.
Definition: residual_based_implicit_time_scheme.h:77
Element::DofsVectorType DofsVectorType
DoF vector type definition.
Definition: residual_based_implicit_time_scheme.h:68
KRATOS_CLASS_POINTER_DEFINITION(ResidualBasedImplicitTimeScheme)
virtual void AddDynamicsToRHS(Condition &rCurrentCondition, LocalSystemVectorType &rRHSContribution, LocalSystemMatrixType &D, LocalSystemMatrixType &M, const ProcessInfo &rCurrentProcessInfo)
It adds the dynamic RHS contribution of the condition RHS = fext - M*an0 - D*vn0 - K*dn0.
Definition: residual_based_implicit_time_scheme.h:406
std::size_t IndexType
Index type definition.
Definition: residual_based_implicit_time_scheme.h:89
BaseType::DofsArrayType DofsArrayType
DoF array type definition.
Definition: residual_based_implicit_time_scheme.h:66
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.
Definition: residual_based_implicit_time_scheme.h:202
BaseType::TSystemVectorType TSystemVectorType
Vector type definition.
Definition: residual_based_implicit_time_scheme.h:75
ModelPart::ElementsContainerType ElementsArrayType
Elements containers definition.
Definition: residual_based_implicit_time_scheme.h:84
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
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
virtual Parameters ValidateAndAssignParameters(Parameters ThisParameters, const Parameters DefaultParameters) const
This method validate and assign default parameters.
Definition: scheme.h:773
virtual void InitializeSolutionStep(ModelPart &rModelPart, TSystemMatrixType &A, TSystemVectorType &Dx, TSystemVectorType &b)
Function called once at the beginning of each solution step.
Definition: scheme.h:272
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
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....
Definition: scheme.h:508
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
std::size_t IndexType
The definition of the index type.
Definition: key_hash.h:35
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
delta_time
Definition: generate_frictional_mortar_condition.py:130
def num_threads
Definition: script.py:75
Definition: residual_based_implicit_time_scheme.h:346
std::vector< Matrix > M
Definition: residual_based_implicit_time_scheme.h:347
std::vector< Matrix > D
First derivative matrix (usually mass matrix)
Definition: residual_based_implicit_time_scheme.h:348