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.
dam_added_mass_condition_process.hpp
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: Lorenzo Gracia
11 //
12 //
13 
14 #if !defined(KRATOS_DAM_ADDED_MASS_CONDITION_PROCESS)
15 #define KRATOS_DAM_ADDED_MASS_CONDITION_PROCESS
16 
17 #include <cmath>
18 
19 // Project includes
20 #include "includes/kratos_flags.h"
22 #include "processes/process.h"
23 
24 // Application includes
26 
27 namespace Kratos
28 {
29 
31 {
32 
33  public:
35 
37 
38  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
39 
42  Parameters &rParameters) : Process(Flags()), mrModelPart(rModelPart)
43  {
45 
46  //only include validation with c++11 since raw_literals do not exist in c++03
47  Parameters default_parameters(R"(
48  {
49  "model_part_name":"PLEASE_CHOOSE_MODEL_PART_NAME",
50  "variable_name": "PLEASE_PRESCRIBE_VARIABLE_NAME",
51  "Modify" : true,
52  "Gravity_Direction" : "Y",
53  "Reservoir_Bottom_Coordinate_in_Gravity_Direction" : 0.0,
54  "Spe_weight" : 0.0,
55  "Water_level" : 0.0,
56  "interval":[
57  0.0,
58  0.0
59  ]
60  } )");
61 
62  // Some values need to be mandatorily prescribed since no meaningful default value exist. For this reason try accessing to them
63  // So that an error is thrown if they don't exist
64  rParameters["Reservoir_Bottom_Coordinate_in_Gravity_Direction"];
65  rParameters["variable_name"];
66  rParameters["model_part_name"];
67 
68  // Now validate agains defaults -- this also ensures no type mismatch
69  rParameters.ValidateAndAssignDefaults(default_parameters);
70 
71  mVariableName = rParameters["variable_name"].GetString();
72  mGravityDirection = rParameters["Gravity_Direction"].GetString();
73  mReferenceCoordinate = rParameters["Reservoir_Bottom_Coordinate_in_Gravity_Direction"].GetDouble();
74  mSpecific = rParameters["Spe_weight"].GetDouble();
75  mWaterLevel = rParameters["Water_level"].GetDouble();
76 
77  KRATOS_CATCH("");
78  }
79 
81 
84 
85  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
86 
87  void ExecuteInitialize() override
88  {
89  KRATOS_TRY;
90 
92  const int nnodes = mrModelPart.GetMesh(0).Nodes().size();
93  int direction;
94  double added_mass;
95 
96  if (mGravityDirection == "X")
97  direction = 0;
98  else if (mGravityDirection == "Y")
99  direction = 1;
100  else
101  direction = 2;
102 
103  if (nnodes != 0)
104  {
105  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.GetMesh(0).NodesBegin();
106 
107 #pragma omp parallel for
108  for (int i = 0; i < nnodes; i++)
109  {
110  ModelPart::NodesContainerType::iterator it = it_begin + i;
111 
112  double y_water = mWaterLevel - (it->Coordinates()[direction]);
113 
114  if (y_water < 0.0)
115  {
116  y_water = 0.0;
117  }
118 
119  added_mass = 0.875 * mSpecific * sqrt(y_water * (mWaterLevel - mReferenceCoordinate));
120 
121  it->FastGetSolutionStepValue(var) = added_mass;
122 
123  if (added_mass > 0.0)
124  {
125  it->FastGetSolutionStepValue(var) = added_mass;
126  }
127  else
128  {
129  it->FastGetSolutionStepValue(var) = 0.0;
130  }
131  }
132  }
133 
134  KRATOS_CATCH("");
135  }
136 
137  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
138 
140  {
141  KRATOS_TRY;
142 
144  const int nnodes = mrModelPart.GetMesh(0).Nodes().size();
145  int direction;
146  double added_mass;
147 
148  if (mGravityDirection == "X")
149  direction = 0;
150  else if (mGravityDirection == "Y")
151  direction = 1;
152  else
153  direction = 2;
154 
155  if (nnodes != 0)
156  {
157  ModelPart::NodesContainerType::iterator it_begin = mrModelPart.GetMesh(0).NodesBegin();
158 
159 #pragma omp parallel for
160  for (int i = 0; i < nnodes; i++)
161  {
162  ModelPart::NodesContainerType::iterator it = it_begin + i;
163 
164  double y_water = mWaterLevel - (it->Coordinates()[direction]);
165 
166  if (y_water < 0.0)
167  {
168  y_water = 0.0;
169  }
170 
171  added_mass = 0.875 * mSpecific * sqrt(y_water * (mWaterLevel - mReferenceCoordinate));
172 
173  it->FastGetSolutionStepValue(var) = added_mass;
174 
175  if (added_mass > 0.0)
176  {
177  it->FastGetSolutionStepValue(var) = added_mass;
178  }
179  else
180  {
181  it->FastGetSolutionStepValue(var) = 0.0;
182  }
183  }
184  }
185 
186  KRATOS_CATCH("");
187  }
188 
190  std::string Info() const override
191  {
192  return "DamAddedMassConditionProcess";
193  }
194 
196  void PrintInfo(std::ostream &rOStream) const override
197  {
198  rOStream << "DamAddedMassConditionProcess";
199  }
200 
202  void PrintData(std::ostream &rOStream) const override
203  {
204  }
205 
207 
208  protected:
210 
212  std::string mVariableName;
213  std::string mGravityDirection;
215  double mSpecific;
216  double mWaterLevel;
217 
218  //----------------------------------------------------------------------------------------------------------------------------------------------------------------------------------
219 
220  private:
223 
224 }; //Class
225 
227 inline std::istream &operator>>(std::istream &rIStream,
229 
231 inline std::ostream &operator<<(std::ostream &rOStream,
232  const DamAddedMassConditionProcess &rThis)
233 {
234  rThis.PrintInfo(rOStream);
235  rOStream << std::endl;
236  rThis.PrintData(rOStream);
237 
238  return rOStream;
239 }
240 
241 } /* namespace Kratos.*/
242 
243 #endif /* KRATOS_DAM_ADDED_MASS_CONDITION_PROCESS defined */
Definition: dam_added_mass_condition_process.hpp:31
double mSpecific
Definition: dam_added_mass_condition_process.hpp:215
KRATOS_CLASS_POINTER_DEFINITION(DamAddedMassConditionProcess)
void ExecuteInitialize() override
This function is designed for being called at the beginning of the computations right after reading t...
Definition: dam_added_mass_condition_process.hpp:87
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: dam_added_mass_condition_process.hpp:196
ModelPart & mrModelPart
Member Variables.
Definition: dam_added_mass_condition_process.hpp:211
std::string Info() const override
Turn back information as a string.
Definition: dam_added_mass_condition_process.hpp:190
virtual ~DamAddedMassConditionProcess()
Destructor.
Definition: dam_added_mass_condition_process.hpp:83
DamAddedMassConditionProcess(ModelPart &rModelPart, Parameters &rParameters)
Constructor.
Definition: dam_added_mass_condition_process.hpp:41
std::string mVariableName
Definition: dam_added_mass_condition_process.hpp:212
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: dam_added_mass_condition_process.hpp:202
double mWaterLevel
Definition: dam_added_mass_condition_process.hpp:216
Table< double, double > TableType
Definition: dam_added_mass_condition_process.hpp:36
std::string mGravityDirection
Definition: dam_added_mass_condition_process.hpp:213
double mReferenceCoordinate
Definition: dam_added_mass_condition_process.hpp:214
void ExecuteInitializeSolutionStep() override
This function will be executed at every time step BEFORE performing the solve phase.
Definition: dam_added_mass_condition_process.hpp:139
Definition: flags.h:58
KratosComponents class encapsulates a lookup table for a family of classes in a generic way.
Definition: kratos_components.h:49
NodesContainerType & Nodes()
Definition: mesh.h:346
NodeIterator NodesBegin()
Definition: mesh.h:326
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType & GetMesh(IndexType ThisIndex=0)
Definition: model_part.h:1791
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
double GetDouble() const
This method returns the double contained in the current Parameter.
Definition: kratos_parameters.cpp:657
void ValidateAndAssignDefaults(const Parameters &rDefaultParameters)
This function is designed to verify that the parameters under testing match the form prescribed by th...
Definition: kratos_parameters.cpp:1306
std::string GetString() const
This method returns the string contained in the current Parameter.
Definition: kratos_parameters.cpp:684
size_type size() const
Returns the number of elements in the container.
Definition: pointer_vector_set.h:502
The base class for all processes in Kratos.
Definition: process.h:49
Definition: table.h:435
#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
int nnodes
Definition: sensitivityMatrix.py:24
integer i
Definition: TensorModule.f:17