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.
calculate_wave_height_utility.h
Go to the documentation of this file.
1 // | / |
2 // ' / __| _` | __| _ \ __|
3 // . \ | ( | | ( |\__ \.
4 // _|\_\_| \__,_|\__|\___/ ____/
5 // Multi-Physics PfemFluidDynamics Application
6 //
7 // License: BSD License
8 // Kratos default license: kratos/license.txt
9 //
10 // Main authors: Alejandro Cornejo Velazquez
11 // Miguel Maso Sotomayor
12 //
13 
14 #ifndef KRATOS_CALCULATE_WAVE_HEIGHT_UTILITY_H_INCLUDED
15 #define KRATOS_CALCULATE_WAVE_HEIGHT_UTILITY_H_INCLUDED
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
22 #include "includes/model_part.h"
24 
25 namespace Kratos
26 {
27 
30 
37 class KRATOS_API(PFEM_FLUID_DYNAMICS_APPLICATION) CalculateWaveHeightUtility
38 {
39 public:
42 
44 
48 
49  typedef Node NodeType;
50 
52 
56 
60  CalculateWaveHeightUtility(ModelPart &rThisModelPart, Parameters ThisParameters);
61 
66 
70 
74  double Calculate(const array_1d<double,3>& rCoordinates) const;
75 
79 
81  std::string Info() const
82  {
83  return "CalculateWaveHeightUtility";
84  }
85 
87  void PrintInfo(std::ostream &rOStream) const
88  {
89  rOStream << "CalculateWaveHeightUtility";
90  }
91 
93  void PrintData(std::ostream &rOStream) const
94  {
95  }
96 
98 
99 protected:
100 
103 
104  // /// Copy constructor.
105  // CalculateWaveHeightUtility(CalculateWaveHeightUtility const &rOther);
106 
110 
111 
115 
116 
118 
119 private:
120 
123 
124  ModelPart &mrModelPart;
125  array_1d<double,3> mDirection;
126  array_1d<double,3> mCoordinates;
127  double mMeanWaterLevel;
128  bool mUseLocalElementSize;
129  bool mUseNearestNode;
130  double mRelativeRadius;
131  double mAbsoluteRadius;
132  double mMeanElementSize;
133 
137 
141  double CalculateAverage(const array_1d<double,3>& rCoordinates, double SearchRadius) const;
142 
143 
147  double CalculateAverage(const array_1d<double,3>& rCoordinates) const;
148 
152  double CalculateNearest(const array_1d<double,3>& rCoordinates) const;
153 
157 
158  // /// Assignment operator.
159  // CalculateWaveHeightUtility &operator=(CalculateWaveHeightUtility const &rOther);
160 
164 
165 
167 
168 }; // Class CalculateWaveHeightUtility
169 
173 
175 inline std::ostream &operator<<(std::ostream &rOStream,
176  const CalculateWaveHeightUtility &rThis)
177 {
178  rThis.PrintInfo(rOStream);
179  rOStream << std::endl;
180  rThis.PrintData(rOStream);
181 
182  return rOStream;
183 }
184 
186 
187 } // namespace Kratos.
188 
189 #endif // KRATOS_CALCULATE_WAVE_HEIGHT_UTILITY_H_INCLUDED defined
This function computes the wave height at a given point.
Definition: calculate_wave_height_utility.h:38
~CalculateWaveHeightUtility()
Destructor.
Definition: calculate_wave_height_utility.h:65
std::string Info() const
Turn back information as a string.
Definition: calculate_wave_height_utility.h:81
Node NodeType
Definition: calculate_wave_height_utility.h:49
void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: calculate_wave_height_utility.h:87
void PrintData(std::ostream &rOStream) const
Print information about this object.
Definition: calculate_wave_height_utility.h:93
Geometry< NodeType > GeometryType
Definition: calculate_wave_height_utility.h:51
KRATOS_CLASS_POINTER_DEFINITION(CalculateWaveHeightUtility)
Geometry base class.
Definition: geometry.h:71
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
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
TDataType Calculate(GeometryType &dummy, const Variable< TDataType > &rVariable)
Definition: add_geometries_to_python.cpp:103
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432