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_slip.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 
14 #ifndef KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_SLIP_H
15 #define KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_SLIP_H
16 
17 /* System includes */
18 
19 
20 /* External includes */
21 
22 
23 /* Project includes */
24 #include "includes/define.h"
25 #include "includes/model_part.h"
27 //#include "includes/variables.h"
31 
32 namespace Kratos
33 {
34 
37 
38 
42 
44 
45 
48 
49 
53 
54 
55 
59 
61 
66 template<class TSparseSpace,
67  class TDenseSpace //= DenseSpace<double>
68  >
70 {
71 
72 public:
75 
77 
79 
81 
83 
84  typedef typename BaseType::TDataType TDataType;
85 
87 
89 
91 
96 
100 
105  {
106  }
107 
113  : BaseType()
114  {
115  // Validate default parameters
116  Parameters default_parameters = Parameters(R"(
117  {
118  "name" : "ResidualBasedIncrementalUpdateStaticSchemeSlip",
119  "domain_size" : 3,
120  "block_size" : 3
121  })" );
122  ThisParameters.ValidateAndAssignDefaults(default_parameters);
123 
124  const int domain_size = ThisParameters["domain_size"].GetInt();
125  const int block_size = ThisParameters["block_size"].GetInt();
126  mpRotationTool = Kratos::make_shared<RotationToolType>(domain_size, block_size, SLIP);
127  }
128 
130 
133  explicit ResidualBasedIncrementalUpdateStaticSchemeSlip(unsigned int DomainSize,
134  unsigned int BlockSize):
135  BaseType(),
136  mpRotationTool(Kratos::make_shared<RotationToolType>(DomainSize,BlockSize,SLIP))
137  {}
138 
140 
144  BaseType(),
145  mpRotationTool(pRotationTool)
146  {}
147 
150 
154 
155 
159 
164  typename BaseSchemeType::Pointer Create(Parameters ThisParameters) const override
165  {
166  return Kratos::make_shared<ClassType>(ThisParameters);
167  }
168 
170  void Update(ModelPart& r_model_part,
171  DofsArrayType& rDofSet,
173  TSystemVectorType& Dx,
174  TSystemVectorType& b) override
175  {
176  KRATOS_TRY;
177 
178  mpRotationTool->RotateVelocities(r_model_part);
179 
180  BaseType::Update(r_model_part,rDofSet,A,Dx,b);
181 
182  mpRotationTool->RecoverVelocities(r_model_part);
183 
184  KRATOS_CATCH("");
185  }
186 
187 
189  void CalculateSystemContributions(Element& rCurrentElement,
190  LocalSystemMatrixType& LHS_Contribution,
191  LocalSystemVectorType& RHS_Contribution,
193  const ProcessInfo& CurrentProcessInfo) override
194  {
195  KRATOS_TRY;
196 
197  BaseType::CalculateSystemContributions(rCurrentElement,LHS_Contribution,RHS_Contribution,EquationId,CurrentProcessInfo);
198 
199  mpRotationTool->Rotate(LHS_Contribution,RHS_Contribution,rCurrentElement.GetGeometry());
200  mpRotationTool->ApplySlipCondition(LHS_Contribution,RHS_Contribution,rCurrentElement.GetGeometry());
201 
202  KRATOS_CATCH("");
203  }
204 
206  void CalculateRHSContribution(Element& rCurrentElement,
207  LocalSystemVectorType& RHS_Contribution,
209  const ProcessInfo& CurrentProcessInfo) override
210  {
211  KRATOS_TRY;
212 
213  BaseType::CalculateRHSContribution(rCurrentElement,RHS_Contribution,EquationId,CurrentProcessInfo);
214 
215  mpRotationTool->Rotate(RHS_Contribution,rCurrentElement.GetGeometry());
216  mpRotationTool->ApplySlipCondition(RHS_Contribution,rCurrentElement.GetGeometry());
217 
218  KRATOS_CATCH("");
219  }
220 
222  void CalculateLHSContribution(Element& rCurrentElement,
223  LocalSystemMatrixType& LHS_Contribution,
225  const ProcessInfo& CurrentProcessInfo) override
226  {
227  KRATOS_TRY;
228 
229  BaseType::CalculateLHSContribution(rCurrentElement,LHS_Contribution,EquationId,CurrentProcessInfo);
230 
231  LocalSystemVectorType Temp = ZeroVector(LHS_Contribution.size1());
232  mpRotationTool->Rotate(LHS_Contribution,Temp,rCurrentElement.GetGeometry());
233  mpRotationTool->ApplySlipCondition(LHS_Contribution,Temp,rCurrentElement.GetGeometry());
234 
235  KRATOS_CATCH("");
236  }
237 
238 
240  void CalculateSystemContributions(Condition& rCurrentCondition,
241  LocalSystemMatrixType& LHS_Contribution,
242  LocalSystemVectorType& RHS_Contribution,
244  const ProcessInfo& CurrentProcessInfo) override
245  {
246  KRATOS_TRY;
247 
248  BaseType::CalculateSystemContributions(rCurrentCondition,LHS_Contribution,RHS_Contribution,EquationId,CurrentProcessInfo);
249 
250  mpRotationTool->Rotate(LHS_Contribution,RHS_Contribution,rCurrentCondition.GetGeometry());
251  mpRotationTool->ApplySlipCondition(LHS_Contribution,RHS_Contribution,rCurrentCondition.GetGeometry());
252 
253  KRATOS_CATCH("");
254  }
255 
257  void CalculateRHSContribution(Condition& rCurrentCondition,
258  LocalSystemVectorType& RHS_Contribution,
260  const ProcessInfo& CurrentProcessInfo) override
261  {
262  KRATOS_TRY;
263 
264  BaseType::CalculateRHSContribution(rCurrentCondition,RHS_Contribution,EquationId,CurrentProcessInfo);
265 
266  mpRotationTool->Rotate(RHS_Contribution,rCurrentCondition.GetGeometry());
267  mpRotationTool->ApplySlipCondition(RHS_Contribution,rCurrentCondition.GetGeometry());
268 
269  KRATOS_CATCH("");
270  }
271 
276  static std::string Name()
277  {
278  return "static_slip_scheme";
279  }
280 
284 
285 
289 
291  std::string Info() const override
292  {
293  return "ResidualBasedIncrementalUpdateStaticSchemeSlip";
294  }
295 
297  void PrintInfo(std::ostream& rOStream) const override
298  {
299  rOStream << Info();
300  }
301 
303  void PrintData(std::ostream& rOStream) const override
304  {
305  rOStream << Info();
306  }
307 
311 
312 
314 
315 protected:
316 
319 
320 
324 
325 
329 
330 
334 
335 
339 
340 
344 
345 
349 
350 
352 
353 private:
354 
357 
358 
362 
364  RotationToolPointerType mpRotationTool;
365 
369 
370 
374 
375 
379 
380 
384 
385 
389 
390 
394 
395 
397 
398 }; // class
399 
403 
405 
407 
408 } // namespace Kratos
409 
410 #endif // KRATOS_RESIDUALBASED_INCREMENTALUPDATE_STATIC_SCHEME_SLIP_H
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