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_process.hpp
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:
9 // kratos/license.txt
10 //
11 // Main authors: Alejandro Cornejo Velazquez
12 //
13 
14 #if !defined(KRATOS_CALCULATE_WAVE_HEIGHT_PROCESS_H_INCLUDED)
15 #define KRATOS_CALCULATE_WAVE_HEIGHT_PROCESS_H_INCLUDED
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
22 #include "includes/model_part.h"
23 #include "processes/process.h"
24 #include <fstream>
25 #include <iostream>
26 
27 namespace Kratos
28 {
29 
32 
34 
37  {
38  public:
41 
44 
49  const int HeightDirection,
50  const int PlaneDirection,
51  const double PlaneCoordinates = 0.0,
52  const double HeightReference = 0.0,
53  const double Tolerance = 1.0e-2,
54  const std::string OutputFileName = "WaveHeight",
55  const double TimeInterval = 0.0) : mrModelPart(rModelPart), mHeightDirection(HeightDirection),
56  mPlaneDirection(PlaneDirection), mPlaneCoordinates(PlaneCoordinates),
57  mHeightReference(HeightReference), mTolerance(Tolerance),
58  mOutputFileName(OutputFileName), mTimeInterval(TimeInterval)
59  {
60  std::ofstream my_file;
61  const std::string file_name = mOutputFileName + ".txt";
62  my_file.open(file_name, std::ios_base::trunc);
63  my_file << " TIME Wave Height" << std::endl;
64  my_file.close();
65  }
66 
69 
73 
75  void operator()()
76  {
77  Execute();
78  }
79 
83 
85  void Execute() override
86  {
88  const double time = mrModelPart.GetProcessInfo()[TIME];
89  const int step = mrModelPart.GetProcessInfo()[STEP];
90 
91  if (time - mPreviousPlotTime > mTimeInterval || step == 1) {
92  // We loop over the nodes...
93  const auto it_node_begin = mrModelPart.NodesBegin();
95  std::vector<double> max_vector(num_threads, -1.0);
96 
97  #pragma omp parallel for
98  for (int i = 0; i < static_cast<int>(mrModelPart.Nodes().size()); i++) {
99  auto it_node = it_node_begin + i;
100 
101  const int thread_id = OpenMPUtils::ThisThread();
102  const auto &r_node_coordinates = it_node->Coordinates();
103  if (it_node->IsNot(ISOLATED) &&
104  it_node->Is(FREE_SURFACE) &&
105  r_node_coordinates(mPlaneDirection) < mPlaneCoordinates + mTolerance &&
106  r_node_coordinates(mPlaneDirection) > mPlaneCoordinates - mTolerance)
107  {
108  const double height = r_node_coordinates(mHeightDirection);
109 
110  const double wave_height = std::abs(height - mHeightReference);
111  if (wave_height > max_vector[thread_id])
112  max_vector[thread_id] = wave_height;
113  }
114  }
115  const double max_height = *std::max_element(max_vector.begin(), max_vector.end());
116 
117  // We open the file where we print the wave height values
118  if (max_height > -1.0) {
119  std::ofstream my_file;
120  const std::string file_name = mOutputFileName + ".txt";
121  my_file.open(file_name, std::ios_base::app);
122  my_file << " " + std::to_string(time) + " " + std::to_string(max_height - mHeightReference) << std::endl;
123  mPreviousPlotTime = time;
124  }
125 
126  }
127  KRATOS_CATCH("");
128  }
129 
133 
137 
141 
143  std::string Info() const override
144  {
145  return "CalculateWaveHeightProcess";
146  }
147 
149  void PrintInfo(std::ostream &rOStream) const override
150  {
151  rOStream << "CalculateWaveHeightProcess";
152  }
153 
158 
159  protected:
168 
171 
185 
186  private:
192 
193  ModelPart &mrModelPart;
194  int mHeightDirection;
195  int mPlaneDirection;
196 
197  double mPlaneCoordinates;
198  double mHeightReference;
199  double mTolerance;
200  std::string mOutputFileName;
201  double mTimeInterval;
202  double mPreviousPlotTime = 0.0;
203 
204 
211 
213  CalculateWaveHeightProcess &operator=(CalculateWaveHeightProcess const &rOther);
214 
224 
225  }; // Class CalculateWaveHeightProcess
226 
228 
231 
235 
237  inline std::istream &operator>>(std::istream &rIStream,
239 
241  inline std::ostream &operator<<(std::ostream &rOStream,
242  const CalculateWaveHeightProcess &rThis)
243  {
244  rThis.PrintInfo(rOStream);
245  rOStream << std::endl;
246  rThis.PrintData(rOStream);
247 
248  return rOStream;
249  }
251 
252 } // namespace Kratos.
253 
254 #endif // KRATOS_TRANSFER_MODEL_PART_ELEMENTS_PROCESS_H_INCLUDED defined
The base class for assigning a value to scalar variables or array_1d components processes in Kratos.
Definition: calculate_wave_height_process.hpp:37
std::string Info() const override
Turn back information as a string.
Definition: calculate_wave_height_process.hpp:143
virtual ~CalculateWaveHeightProcess()
Destructor.
Definition: calculate_wave_height_process.hpp:68
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: calculate_wave_height_process.hpp:149
CalculateWaveHeightProcess(ModelPart &rModelPart, const int HeightDirection, const int PlaneDirection, const double PlaneCoordinates=0.0, const double HeightReference=0.0, const double Tolerance=1.0e-2, const std::string OutputFileName="WaveHeight", const double TimeInterval=0.0)
Definition: calculate_wave_height_process.hpp:48
KRATOS_CLASS_POINTER_DEFINITION(CalculateWaveHeightProcess)
Pointer definition of CalculateWaveHeightProcess.
CalculateWaveHeightProcess(CalculateWaveHeightProcess const &rOther)
Copy constructor.
void Execute() override
Execute method is used to execute the CalculateWaveHeightProcess.
Definition: calculate_wave_height_process.hpp:85
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: calculate_wave_height_process.hpp:75
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ProcessInfo & GetProcessInfo()
Definition: model_part.h:1746
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
static int ThisThread()
Wrapper for omp_get_thread_num().
Definition: openmp_utils.h:108
static int GetNumThreads()
Returns the current number of threads.
Definition: parallel_utilities.cpp:34
The base class for all processes in Kratos.
Definition: process.h:49
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: process.h:210
#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
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
time
Definition: face_heat.py:85
int step
Definition: face_heat.py:88
def num_threads
Definition: script.py:75
string file_name
Definition: sp_statistics.py:6
integer i
Definition: TensorModule.f:17