14 #if !defined(KRATOS_CALCULATE_WAVE_HEIGHT_PROCESS_H_INCLUDED)
15 #define KRATOS_CALCULATE_WAVE_HEIGHT_PROCESS_H_INCLUDED
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)
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;
91 if (
time - mPreviousPlotTime > mTimeInterval ||
step == 1) {
93 const auto it_node_begin = mrModelPart.
NodesBegin();
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;
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)
108 const double height = r_node_coordinates(mHeightDirection);
110 const double wave_height = std::abs(height - mHeightReference);
111 if (wave_height > max_vector[thread_id])
112 max_vector[thread_id] = wave_height;
115 const double max_height = *std::max_element(max_vector.begin(), max_vector.end());
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;
143 std::string
Info()
const override
145 return "CalculateWaveHeightProcess";
151 rOStream <<
"CalculateWaveHeightProcess";
194 int mHeightDirection;
197 double mPlaneCoordinates;
198 double mHeightReference;
200 std::string mOutputFileName;
201 double mTimeInterval;
202 double mPreviousPlotTime = 0.0;
245 rOStream << std::endl;
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