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.
residualbased_incrementalupdate_static_scheme.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ `
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Riccardo Rossi
11 //
12 
13 #if !defined(KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_H )
14 #define KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_H
15 
16 /* System includes */
17 
18 /* External includes */
19 
20 /* Project includes */
22 #include "includes/variables.h"
24 
25 namespace Kratos
26 {
41 
52 template<class TSparseSpace,
53  class TDenseSpace //= DenseSpace<double>
54  >
56  : public Scheme<TSparseSpace,TDenseSpace>
57 {
58 
59 public:
62 
65 
68 
69  // The current class definition
71 
74 
76  typedef typename BaseType::TDataType TDataType;
85 
90 
93 
97 
103  : BaseType()
104  {
105  // Validate and assign defaults
106  ThisParameters = this->ValidateAndAssignParameters(ThisParameters, this->GetDefaultParameters());
107  this->AssignSettings(ThisParameters);
108  }
109 
113  : BaseType()
114  {}
115 
119  :BaseType(rOther)
120  {
121  }
122 
126 
130 
134 
139  typename BaseType::Pointer Create(Parameters ThisParameters) const override
140  {
141  return Kratos::make_shared<ClassType>(ThisParameters);
142  }
143 
152  void Update(
153  ModelPart& rModelPart,
154  DofsArrayType& rDofSet,
155  TSystemMatrixType& rA,
156  TSystemVectorType& rDx,
158  ) override
159  {
160  KRATOS_TRY
161 
162  mpDofUpdater->UpdateDofs(rDofSet, rDx);
163 
164  KRATOS_CATCH("")
165  }
166 
174  void Predict(
175  ModelPart& rModelPart,
176  DofsArrayType& rDofSet,
177  TSystemMatrixType& rA,
178  TSystemVectorType& rDx,
180  ) override
181  {
182  KRATOS_TRY
183 
184  KRATOS_CATCH("")
185  }
186 
197  Element& rCurrentElement,
198  LocalSystemMatrixType& rLHSContribution,
199  LocalSystemVectorType& rRHSContribution,
200  EquationIdVectorType& rEquationId,
201  const ProcessInfo& rCurrentProcessInfo
202  ) override
203  {
204  KRATOS_TRY
205 
206  rCurrentElement.CalculateLocalSystem(rLHSContribution,rRHSContribution, rCurrentProcessInfo);
207 
208  rCurrentElement.EquationIdVector(rEquationId, rCurrentProcessInfo);
209 
210  KRATOS_CATCH("")
211  }
212 
222  Condition& rCurrentCondition,
223  LocalSystemMatrixType& rLHSContribution,
224  LocalSystemVectorType& rRHSContribution,
225  EquationIdVectorType& rEquationId,
226  const ProcessInfo& rCurrentProcessInfo
227  ) override
228  {
229  KRATOS_TRY
230 
231  rCurrentCondition.CalculateLocalSystem(rLHSContribution, rRHSContribution, rCurrentProcessInfo);
232 
233  rCurrentCondition.EquationIdVector(rEquationId, rCurrentProcessInfo);
234 
235  KRATOS_CATCH("")
236  }
237 
246  Element& rCurrentElement,
247  LocalSystemVectorType& rRHSContribution,
248  EquationIdVectorType& rEquationId,
249  const ProcessInfo& rCurrentProcessInfo
250  ) override
251  {
252  KRATOS_TRY
253 
254  rCurrentElement.CalculateRightHandSide(rRHSContribution, rCurrentProcessInfo);
255 
256  rCurrentElement.EquationIdVector(rEquationId, rCurrentProcessInfo);
257 
258  KRATOS_CATCH("")
259  }
260 
269  Condition& rCurrentCondition,
270  LocalSystemVectorType& rRHSContribution,
271  EquationIdVectorType& rEquationId,
272  const ProcessInfo& rCurrentProcessInfo
273  ) override
274  {
275  KRATOS_TRY
276 
277  rCurrentCondition.CalculateRightHandSide(rRHSContribution, rCurrentProcessInfo);
278 
279  rCurrentCondition.EquationIdVector(rEquationId, rCurrentProcessInfo);
280 
281  KRATOS_CATCH("")
282  }
283 
292  Element& rCurrentElement,
293  LocalSystemMatrixType& rLHSContribution,
294  EquationIdVectorType& rEquationId,
295  const ProcessInfo& rCurrentProcessInfo
296  ) override
297  {
298  KRATOS_TRY
299 
300  rCurrentElement.CalculateLeftHandSide(rLHSContribution, rCurrentProcessInfo);
301 
302  rCurrentElement.EquationIdVector(rEquationId, rCurrentProcessInfo);
303 
304  KRATOS_CATCH("")
305  }
306 
310  void Clear() override
311  {
312  this->mpDofUpdater->Clear();
313  }
314 
320  {
321  Parameters default_parameters = Parameters(R"(
322  {
323  "name" : "static_scheme"
324  })");
325 
326  // Getting base class default parameters
327  const Parameters base_default_parameters = BaseType::GetDefaultParameters();
328  default_parameters.RecursivelyAddMissingParameters(base_default_parameters);
329  return default_parameters;
330  }
331 
336  static std::string Name()
337  {
338  return "static_scheme";
339  }
340 
344 
348 
352 
354  std::string Info() const override
355  {
356  return "ResidualBasedIncrementalUpdateStaticScheme";
357  }
358 
360  void PrintInfo(std::ostream& rOStream) const override
361  {
362  rOStream << Info();
363  }
364 
366  void PrintData(std::ostream& rOStream) const override
367  {
368  rOStream << Info();
369  }
370 
374 
376 
377 protected:
380 
384 
388 
392 
396 
400 
404 
406 
407 private:
410 
414 
415  typename TSparseSpace::DofUpdaterPointerType mpDofUpdater = TSparseSpace::CreateDofUpdater();
416 
420 
424 
428 
432 
436 
438 
439 }; // Class ResidualBasedIncrementalUpdateStaticScheme
440 } // namespace Kratos
441 
442 #endif /* KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_H defined */
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