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.
wave_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 // Main authors: Miguel Maso Sotomayor
11 //
12 
13 #ifndef KRATOS_WAVE_CONDITION_H_INCLUDED
14 #define KRATOS_WAVE_CONDITION_H_INCLUDED
15 
16 // System includes
17 
18 
19 // External includes
20 
21 
22 // Project includes
23 #include "includes/define.h"
24 #include "includes/condition.h"
25 #include "includes/serializer.h"
26 
27 namespace Kratos
28 {
31 
34 
38 
42 
49 template<std::size_t TNumNodes>
50 class WaveCondition : public Condition
51 {
52 public:
55 
56  typedef std::size_t IndexType;
57 
58  typedef Node NodeType;
59 
61 
63 
65 
68 
72 
76 
81 
85  WaveCondition(IndexType NewId, const NodesArrayType& ThisNodes) : Condition(NewId, ThisNodes){}
86 
90  WaveCondition(IndexType NewId, GeometryType::Pointer pGeometry) : Condition(NewId, pGeometry){}
91 
95  WaveCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties) : Condition(NewId, pGeometry, pProperties){}
96 
100  ~ WaveCondition() override {};
101 
105 
113  Condition::Pointer Create(IndexType NewId, NodesArrayType const& ThisNodes, PropertiesType::Pointer pProperties) const override
114  {
115  return Kratos::make_intrusive<WaveCondition<TNumNodes>>(NewId, this->GetGeometry().Create(ThisNodes), pProperties);
116  }
117 
125  Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
126  {
127  return Kratos::make_intrusive<WaveCondition<TNumNodes>>(NewId, pGeom, pProperties);
128  }
129 
136  Condition::Pointer Clone(IndexType NewId, NodesArrayType const& ThisNodes) const override
137  {
138  Condition::Pointer p_new_elem = Create(NewId, this->GetGeometry().Create(ThisNodes), this->pGetProperties());
139  p_new_elem->SetData(this->GetData());
140  p_new_elem->Set(Flags(*this));
141  return p_new_elem;
142  }
143 
148  int Check(const ProcessInfo& rCurrentProcessInfo) const override;
149 
155  void EquationIdVector(EquationIdVectorType& rResult, const ProcessInfo& rCurrentProcessInfo) const override;
156 
162  void GetDofList(DofsVectorType& rConditionalDofList, const ProcessInfo& rCurrentProcessInfo) const override;
163 
167  void GetValuesVector(Vector& rValues, int Step = 0) const override;
168 
172  void GetFirstDerivativesVector(Vector& rValues, int Step = 0) const override;
173 
177  void GetSecondDerivativesVector(Vector& rValues, int Step = 0) const override;
178 
185  void Calculate(
186  const Variable<array_1d<double,3>>& rVariable,
187  array_1d<double,3>& rOutput,
188  const ProcessInfo& rCurrentProcessInfo) override;
189 
195  void CalculateRightHandSide(VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
196 
203  void CalculateLocalSystem(MatrixType& rLeftHandSideMatrix, VectorType& rRightHandSideVector, const ProcessInfo& rCurrentProcessInfo) override;
204 
210  void CalculateMassMatrix(MatrixType& rMassMatrix, const ProcessInfo& rCurrentProcessInfo) override;
211 
215 
220  const Parameters GetSpecifications() const override;
221 
225 
229  std::string Info() const override
230  {
231  return "WaveCondition";
232  }
233 
237  void PrintInfo(std::ostream& rOStream) const override
238  {
239  rOStream << Info() << " : " << Id();
240  }
241 
245  void PrintData(std::ostream& rOStream) const override
246  {
247  rOStream << GetGeometry();
248  }
249 
253 
254 
256 
257 protected:
260 
261  static constexpr IndexType mLocalSize = 3 * TNumNodes;
262 
266 
268  {
270  double stab_factor;
272  double length;
273  double gravity;
274 
275  double depth;
276  double height;
278 
279  double v_neumann;
280  double h_dirichlet;
282 
288 
290  };
291 
295 
296  virtual const Variable<double>& GetUnknownComponent(int Index) const;
297 
299 
300  virtual void InitializeData(
301  ConditionData& rData,
302  const ProcessInfo& rProcessInfo);
303 
304  virtual void CalculateGaussPointData(
305  ConditionData& rData,
306  const IndexType PointIndex,
307  const array_1d<double,TNumNodes>& rN);
308 
309  static void CalculateGeometryData(
310  const GeometryType& rGeometry,
311  Vector &rGaussWeights,
312  Matrix &rNContainer);
313 
314  static const array_1d<double,3> VectorProduct(
315  const array_1d<array_1d<double,3>,TNumNodes>& rV,
316  const array_1d<double,TNumNodes>& rN);
317 
318 
319  void AddFluxTerms(
320  LocalVectorType& rVector,
321  const ConditionData& rData,
322  const array_1d<double,TNumNodes>& rN,
323  const double Weight = 1.0);
324 
325  void AddMassTerms(
326  LocalMatrixType& rMatrix,
327  const ConditionData& rData,
328  const array_1d<double,TNumNodes>& rN,
329  const double Weight);
330 
332 
333 private:
336 
337 
341 
342 
346 
347  friend class Serializer;
348 
349  void save(Serializer& rSerializer) const override
350  {
352  }
353 
354  void load(Serializer& rSerializer) override
355  {
357  }
358 
362 
363 
365 
366 }; // Class WaveCondition
367 
371 
372 
376 
377 
379 
381 
382 } // namespace Kratos.
383 
384 #endif // KRATOS_WAVE_CONDITION_H_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
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
PropertiesType::Pointer pGetProperties()
returns the pointer to the property of the condition. Does not throw an error, to allow copying of co...
Definition: condition.h:964
std::size_t IndexType
Definition: flags.h:74
Flags()
Default constructor.
Definition: flags.h:119
GeometryType & GetGeometry()
Returns the reference of the geometry.
Definition: geometrical_object.h:158
DataValueContainer & GetData()
Definition: geometrical_object.h:212
Geometry base class.
Definition: geometry.h:71
virtual Pointer Create(PointsArrayType const &rThisPoints) const
Creates a new geometry pointer.
Definition: geometry.h:813
IndexType Id() const
Definition: indexed_object.h:107
Definition: amatrix_interface.h:41
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
Implementation of a condition for shallow water waves problems.
Definition: wave_condition.h:51
virtual const Variable< double > & GetUnknownComponent(int Index) const
Definition: wave_condition.cpp:194
WaveCondition(IndexType NewId, GeometryType::Pointer pGeometry, PropertiesType::Pointer pProperties)
Constructor using Geometry and Properties.
Definition: wave_condition.h:95
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: wave_condition.h:245
void Calculate(const Variable< array_1d< double, 3 >> &rVariable, array_1d< double, 3 > &rOutput, const ProcessInfo &rCurrentProcessInfo) override
Access for variables on Integration points.
Definition: wave_condition.cpp:159
int Check(const ProcessInfo &rCurrentProcessInfo) const override
Check that all required data containers are properly initialized and registered in Kratos.
Definition: wave_condition.cpp:48
void AddFluxTerms(LocalVectorType &rVector, const ConditionData &rData, const array_1d< double, TNumNodes > &rN, const double Weight=1.0)
Definition: wave_condition.cpp:331
WaveCondition()
Default constructor.
Definition: wave_condition.h:80
void CalculateLocalSystem(MatrixType &rLeftHandSideMatrix, VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Calculate the conditional contribution to the problem.
Definition: wave_condition.cpp:379
static constexpr IndexType mLocalSize
Definition: wave_condition.h:261
array_1d< double, 3 *TNumNodes > LocalVectorType
Definition: wave_condition.h:60
WaveCondition(IndexType NewId, const NodesArrayType &ThisNodes)
Constructor using an array of nodes.
Definition: wave_condition.h:85
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: wave_condition.h:237
~ WaveCondition() override
Destructor.
Definition: wave_condition.h:100
void GetDofList(DofsVectorType &rConditionalDofList, const ProcessInfo &rCurrentProcessInfo) const override
Fill given array with containing the condition's degrees of freedom.
Definition: wave_condition.cpp:97
std::string Info() const override
Turn back information as a string.
Definition: wave_condition.h:229
Node NodeType
Definition: wave_condition.h:58
KRATOS_CLASS_INTRUSIVE_POINTER_DEFINITION(WaveCondition)
Pointer definition.
virtual void InitializeData(ConditionData &rData, const ProcessInfo &rProcessInfo)
Definition: wave_condition.cpp:257
static void CalculateGeometryData(const GeometryType &rGeometry, Vector &rGaussWeights, Matrix &rNContainer)
Definition: wave_condition.cpp:220
BoundedMatrix< double, 3 *TNumNodes, 3 *TNumNodes > LocalMatrixType
Definition: wave_condition.h:62
WaveCondition(IndexType NewId, GeometryType::Pointer pGeometry)
Constructor using Geometry.
Definition: wave_condition.h:90
void EquationIdVector(EquationIdVectorType &rResult, const ProcessInfo &rCurrentProcessInfo) const override
Fill given vector with the linear system row index for the condition's degrees of freedom.
Definition: wave_condition.cpp:76
const Parameters GetSpecifications() const override
This method provides the specifications/requirements of the element.
Definition: wave_condition.cpp:35
Condition::Pointer Create(IndexType NewId, GeometryType::Pointer pGeom, PropertiesType::Pointer pProperties) const override
Create a new condition pointer.
Definition: wave_condition.h:125
void AddMassTerms(LocalMatrixType &rMatrix, const ConditionData &rData, const array_1d< double, TNumNodes > &rN, const double Weight)
Definition: wave_condition.cpp:357
void GetSecondDerivativesVector(Vector &rValues, int Step=0) const override
Get the second time derivative of variable which defines the degrees of freedom.
Definition: wave_condition.cpp:152
void GetValuesVector(Vector &rValues, int Step=0) const override
Get the variable which defines the degrees of freedom.
Definition: wave_condition.cpp:118
virtual void CalculateGaussPointData(ConditionData &rData, const IndexType PointIndex, const array_1d< double, TNumNodes > &rN)
Definition: wave_condition.cpp:281
static const array_1d< double, 3 > VectorProduct(const array_1d< array_1d< double, 3 >, TNumNodes > &rV, const array_1d< double, TNumNodes > &rN)
Definition: wave_condition.cpp:243
void CalculateMassMatrix(MatrixType &rMassMatrix, const ProcessInfo &rCurrentProcessInfo) override
Calculate the condition mass matrix.
Definition: wave_condition.cpp:413
void GetFirstDerivativesVector(Vector &rValues, int Step=0) const override
Get the time derivative of variable which defines the degrees of freedom.
Definition: wave_condition.cpp:135
void CalculateRightHandSide(VectorType &rRightHandSideVector, const ProcessInfo &rCurrentProcessInfo) override
Calculate the conditional contribution to the problem.
Definition: wave_condition.cpp:367
std::size_t IndexType
Definition: wave_condition.h:56
Condition::Pointer Clone(IndexType NewId, NodesArrayType const &ThisNodes) const override
Create a new condition pointer and clone the previous condition data.
Definition: wave_condition.h:136
GeometryType::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: wave_condition.h:64
virtual LocalVectorType GetUnknownVector(ConditionData &rData)
Definition: wave_condition.cpp:206
Condition::Pointer Create(IndexType NewId, NodesArrayType const &ThisNodes, PropertiesType::Pointer pProperties) const override
Create a new condition pointer.
Definition: wave_condition.h:113
Short class definition.
Definition: array_1d.h:61
#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
def Index()
Definition: hdf5_io_tools.py:38
def load(f)
Definition: ode_solve.py:307
Definition: wave_condition.h:268
double length
Definition: wave_condition.h:272
double relative_dry_height
Definition: wave_condition.h:271
array_1d< double, TNumNodes > nodal_h
Definition: wave_condition.h:284
double v_neumann
Definition: wave_condition.h:279
array_1d< double, TNumNodes > nodal_f
Definition: wave_condition.h:283
double stab_factor
Definition: wave_condition.h:270
double depth
Definition: wave_condition.h:275
array_1d< array_1d< double, 3 >, TNumNodes > nodal_q
Definition: wave_condition.h:287
array_1d< double, 3 > flux
Definition: wave_condition.h:281
array_1d< double, 3 > normal
Definition: wave_condition.h:289
array_1d< array_1d< double, 3 >, TNumNodes > nodal_v
Definition: wave_condition.h:286
double h_dirichlet
Definition: wave_condition.h:280
double height
Definition: wave_condition.h:276
bool integrate_by_parts
Definition: wave_condition.h:269
array_1d< double, TNumNodes > nodal_z
Definition: wave_condition.h:285
array_1d< double, 3 > velocity
Definition: wave_condition.h:277
double gravity
Definition: wave_condition.h:273