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.
support_lagrange_condition.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 
11 #if !defined(KRATOS_SUPPORT_LAGRANGE_CONDITION_H_INCLUDED )
12 #define KRATOS_SUPPORT_LAGRANGE_CONDITION_H_INCLUDED
13 
14 // System includes
15 #include "includes/define.h"
16 #include "includes/condition.h"
17 
18 // External includes
19 
20 // Project includes
22 
23 
24 namespace Kratos
25 {
26 
28 
35  : public Condition
36  {
37  public:
40 
43 
45  typedef std::size_t SizeType;
46  typedef std::size_t IndexType;
47 
51 
54  IndexType NewId,
55  GeometryType::Pointer pGeometry)
56  : Condition(NewId, pGeometry)
57  {};
58 
61  IndexType NewId,
62  GeometryType::Pointer pGeometry,
63  PropertiesType::Pointer pProperties)
64  : Condition(NewId, pGeometry, pProperties)
65  {};
66 
69  : Condition()
70  {};
71 
73  virtual ~SupportLagrangeCondition() = default;
74 
78 
80  Condition::Pointer Create(
81  IndexType NewId,
82  GeometryType::Pointer pGeom,
83  PropertiesType::Pointer pProperties
84  ) const override
85  {
86  return Kratos::make_intrusive<SupportLagrangeCondition>(
87  NewId, pGeom, pProperties);
88  };
89 
91  Condition::Pointer Create(
92  IndexType NewId,
93  NodesArrayType const& ThisNodes,
94  PropertiesType::Pointer pProperties
95  ) const override
96  {
97  return Kratos::make_intrusive< SupportLagrangeCondition >(
98  NewId, GetGeometry().Create(ThisNodes), pProperties);
99  };
100 
104 
107  VectorType& rRightHandSideVector,
108  const ProcessInfo& rCurrentProcessInfo) override
109  {
110  const SizeType mat_size = GetNumberOfNonZeroNodes() * 6;
111 
112  if (rRightHandSideVector.size() != mat_size)
113  rRightHandSideVector.resize(mat_size);
114  noalias(rRightHandSideVector) = ZeroVector(mat_size);
115 
116  MatrixType left_hand_side_matrix;
117 
118  CalculateAll(left_hand_side_matrix, rRightHandSideVector,
119  rCurrentProcessInfo, false, true);
120  }
121 
124  MatrixType& rLeftHandSideMatrix,
125  const ProcessInfo& rCurrentProcessInfo) override
126  {
127  const SizeType mat_size = GetNumberOfNonZeroNodes() * 6;
128 
129  VectorType right_hand_side_vector;
130 
131  if (rLeftHandSideMatrix.size1() != mat_size && rLeftHandSideMatrix.size2())
132  rLeftHandSideMatrix.resize(mat_size, mat_size);
133  noalias(rLeftHandSideMatrix) = ZeroMatrix(mat_size, mat_size);
134 
135  CalculateAll(rLeftHandSideMatrix, right_hand_side_vector,
136  rCurrentProcessInfo, true, false);
137  }
138 
141  MatrixType& rLeftHandSideMatrix,
142  VectorType& rRightHandSideVector,
143  const ProcessInfo& rCurrentProcessInfo) override
144  {
145  const SizeType mat_size = GetNumberOfNonZeroNodes() * 6;
146 
147  if (rRightHandSideVector.size() != mat_size)
148  rRightHandSideVector.resize(mat_size);
149  noalias(rRightHandSideVector) = ZeroVector(mat_size);
150 
151  if (rLeftHandSideMatrix.size1() != mat_size)
152  rLeftHandSideMatrix.resize(mat_size, mat_size);
153  noalias(rLeftHandSideMatrix) = ZeroMatrix(mat_size, mat_size);
154 
155  CalculateAll(rLeftHandSideMatrix, rRightHandSideVector,
156  rCurrentProcessInfo, true, true);
157  }
158 
164  void EquationIdVector(
165  EquationIdVectorType& rResult,
166  const ProcessInfo& rCurrentProcessInfo
167  ) const override;
168 
174  void GetDofList(
175  DofsVectorType& rElementalDofList,
176  const ProcessInfo& rCurrentProcessInfo
177  ) const override;
178 
187  void CalculateAll(
188  MatrixType& rLeftHandSideMatrix,
189  VectorType& rRightHandSideVector,
190  const ProcessInfo& rCurrentProcessInfo,
191  const bool CalculateStiffnessMatrixFlag,
192  const bool CalculateResidualVectorFlag
193  );
194 
195  /* @brief This function computes the determinants of jacobian for
196  * the initial configuration. Required for conservative integration.
197  * @param rGeometry: corresponding geometry
198  * @param rDeterminantOfJacobian: output determinants
199  */
201  const GeometryType& rGeometry,
202  Vector& rDeterminantOfJacobian);
203 
207 
209  std::string Info() const override
210  {
211  std::stringstream buffer;
212  buffer << "\"SupportLagrangeCondition\" #" << Id();
213  return buffer.str();
214  }
215 
217  void PrintInfo(std::ostream& rOStream) const override
218  {
219  rOStream << "\"SupportLagrangeCondition\" #" << Id();
220  }
221 
223  void PrintData(std::ostream& rOStream) const override {
224  pGetGeometry()->PrintData(rOStream);
225  }
226 
228 
229  private:
232 
233  const double shape_function_tolerance = 1e-6;
234 
238 
239  /* @brief checks all shape functions and
240  * returns the number of non zero nodes.
241  */
242  SizeType GetNumberOfNonZeroNodes() const;
243 
247 
248  friend class Serializer;
249 
250  virtual void save(Serializer& rSerializer) const override
251  {
253  }
254 
255  virtual void load(Serializer& rSerializer) override
256  {
258  }
259 
261 
262  }; // Class SupportLagrangeCondition
263 
264 } // namespace Kratos.
265 
266 #endif // KRATOS_SUPPORT_LAGRANGE_CONDITION_H_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
std::size_t SizeType
Definition: condition.h:94
std::vector< std::size_t > EquationIdVectorType
Definition: condition.h:98
Matrix MatrixType
Definition: condition.h:90
std::vector< DofType::Pointer > DofsVectorType
Definition: condition.h:100
std::size_t IndexType
Definition: flags.h:74
GeometryType::Pointer pGetGeometry()
Returns the pointer to the geometry.
Definition: geometrical_object.h:140
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
Geometry base class.
Definition: geometry.h:71
IndexType Id() const
Definition: indexed_object.h:107
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
ProcessInfo holds the current value of different solution parameters.
Definition: process_info.h:59
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Lagrange Multiplier based support condition.
Definition: support_lagrange_condition.h:36
std::string Info() const override
Turn back information as a string.
Definition: support_lagrange_condition.h:209
void DeterminantOfJacobianInitial(const GeometryType &rGeometry, Vector &rDeterminantOfJacobian)
Definition: support_lagrange_condition.cpp:142
void CalculateAll(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo, const bool CalculateStiffnessMatrixFlag, const bool CalculateResidualVectorFlag)
Definition: support_lagrange_condition.cpp:23
std::size_t IndexType
Definition: support_lagrange_condition.h:46
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Sets on rResult the ID's of the element degrees of freedom.
Definition: support_lagrange_condition.cpp:181
SupportLagrangeCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor with Id and geometry.
Definition: support_lagrange_condition.h:53
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create with Id, pointer to geometry and pointer to property.
Definition: support_lagrange_condition.h:91
std::size_t SizeType
Size types.
Definition: support_lagrange_condition.h:45
KRATOS_CLASS_POINTER_DEFINITION(SupportLagrangeCondition)
Counted pointer of SupportLagrangeCondition.
SupportLagrangeCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor with Id, geometry and property.
Definition: support_lagrange_condition.h:60
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Calculates RightHandSide F-vector only.
Definition: support_lagrange_condition.h:106
virtual ~SupportLagrangeCondition()=default
Destructor.
void GetDofList(DofsVectorType &rElementalDofList, const ProcessInfo &rCurrentProcessInfo) const override
Sets on rElementalDofList the degrees of freedom of the considered element geometry.
Definition: support_lagrange_condition.cpp:222
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Calculates LeftHandSide K-Matrix and f-vector.
Definition: support_lagrange_condition.h:140
SupportLagrangeCondition()
Default constructor.
Definition: support_lagrange_condition.h:68
Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
Create with Id, pointer to geometry and pointer to property.
Definition: support_lagrange_condition.h:80
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: support_lagrange_condition.h:223
void CalculateLeftHandSide(MatrixType &rLeftHandSideMatrix, const ProcessInfo &rCurrentProcessInfo) override
Calculates LeftHandSide K-Matrix only.
Definition: support_lagrange_condition.h:123
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: support_lagrange_condition.h:217
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
KratosZeroVector< double > ZeroVector
Definition: amatrix_interface.h:561
KratosZeroMatrix< double > ZeroMatrix
Definition: amatrix_interface.h:559
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
def load(f)
Definition: ode_solve.py:307
e
Definition: run_cpp_mpi_tests.py:31