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.
set_material_properties_to_fluid_nodes_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosPfemFluidApplication $
3 // Created by: $Author: AFranci $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: October 2018 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_SET_MATERIAL_PROPERTIES_TO_FLUID_NODES_PROCESS_H_INCLUDED)
11 #define KRATOS_SET_MATERIAL_PROPERTIES_TO_FLUID_NODES_PROCESS_H_INCLUDED
12 
13 // System includes
14 
15 // External includes
16 
17 // Project includes
18 
20 
23 #include "includes/model_part.h"
24 #include "utilities/openmp_utils.h"
25 #include "utilities/math_utils.h"
27 
29 // Data:
30 // StepData:
31 // Flags: (checked)
32 // (set)
33 // (modified)
34 // (reset)
35 
36 namespace Kratos
37 {
38 
41 
48  typedef std::size_t SizeType;
49 
53 
57 
61 
63 
66  : public MesherProcess
67  {
68  public:
71 
74 
78 
81  : mrModelPart(rModelPart)
82  {
83  }
84 
87  {
88  }
89 
90  void operator()()
91  {
92  Execute();
93  }
94 
98 
99  void Execute() override{
100  KRATOS_TRY
101 
102 #pragma omp parallel
103  {
104 
105  ModelPart::ElementIterator ElemBegin;
108  for (ModelPart::ElementIterator itElem = ElemBegin; itElem != ElemEnd; ++itElem)
109  {
110  ModelPart::PropertiesType &elemProperties = itElem->GetProperties();
111 
112  double flow_index = 1;
113  double yield_shear = 0;
114  double adaptive_exponent = 0;
115  double static_friction = 0;
116  double dynamic_friction = 0;
117  double inertial_number_zero = 0;
118  double grain_diameter = 0;
119  double grain_density = 0;
120  double regularization_coefficient = 0;
121  double friction_angle = 0;
122  double cohesion = 0;
123 
124  double density = elemProperties[DENSITY];
125  double bulk_modulus = elemProperties[BULK_MODULUS];
126  double viscosity = elemProperties[DYNAMIC_VISCOSITY];
127  unsigned int elem_property_id = elemProperties.Id();
128 
129  if (elemProperties.Has(YIELD_SHEAR)) // Bingham model
130  {
131  flow_index = elemProperties[FLOW_INDEX];
132  yield_shear = elemProperties[YIELD_SHEAR];
133  adaptive_exponent = elemProperties[ADAPTIVE_EXPONENT];
134  }
135  else if (elemProperties.Has(INTERNAL_FRICTION_ANGLE)) // Frictional Viscoplastic model
136  {
137  friction_angle = elemProperties[INTERNAL_FRICTION_ANGLE];
138  cohesion = elemProperties[COHESION];
139  adaptive_exponent = elemProperties[ADAPTIVE_EXPONENT];
140  }
141  else if (elemProperties.Has(STATIC_FRICTION)) // Mu(I)-rheology
142  {
143  static_friction = elemProperties[STATIC_FRICTION];
144  dynamic_friction = elemProperties[DYNAMIC_FRICTION];
145  inertial_number_zero = elemProperties[INERTIAL_NUMBER_ZERO];
146  grain_diameter = elemProperties[GRAIN_DIAMETER];
147  grain_density = elemProperties[GRAIN_DENSITY];
148  regularization_coefficient = elemProperties[REGULARIZATION_COEFFICIENT];
149  }
150 
151  Geometry<Node> &rGeom = itElem->GetGeometry();
152  const SizeType NumNodes = rGeom.PointsNumber();
153  for (SizeType i = 0; i < NumNodes; ++i)
154  {
155 
157  {
158  rGeom[i].FastGetSolutionStepValue(PROPERTY_ID) = elem_property_id;
159  }
160 
162  rGeom[i].FastGetSolutionStepValue(BULK_MODULUS) = bulk_modulus;
163 
165  rGeom[i].FastGetSolutionStepValue(DENSITY) = density;
166 
167  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(DYNAMIC_VISCOSITY))
168  rGeom[i].FastGetSolutionStepValue(DYNAMIC_VISCOSITY) = viscosity;
169 
171  rGeom[i].FastGetSolutionStepValue(YIELD_SHEAR) = yield_shear;
172 
174  rGeom[i].FastGetSolutionStepValue(FLOW_INDEX) = flow_index;
175 
176  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(ADAPTIVE_EXPONENT))
177  rGeom[i].FastGetSolutionStepValue(ADAPTIVE_EXPONENT) = adaptive_exponent;
178 
179  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(INTERNAL_FRICTION_ANGLE))
180  rGeom[i].FastGetSolutionStepValue(INTERNAL_FRICTION_ANGLE) = friction_angle;
181 
183  rGeom[i].FastGetSolutionStepValue(COHESION) = cohesion;
184 
185  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(ADAPTIVE_EXPONENT))
186  rGeom[i].FastGetSolutionStepValue(ADAPTIVE_EXPONENT) = adaptive_exponent;
187 
188  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(STATIC_FRICTION))
189  rGeom[i].FastGetSolutionStepValue(STATIC_FRICTION) = static_friction;
190 
191  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(DYNAMIC_FRICTION))
192  rGeom[i].FastGetSolutionStepValue(DYNAMIC_FRICTION) = dynamic_friction;
193 
194  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(INERTIAL_NUMBER_ZERO))
195  rGeom[i].FastGetSolutionStepValue(INERTIAL_NUMBER_ZERO) = inertial_number_zero;
196 
197  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(GRAIN_DIAMETER))
198  rGeom[i].FastGetSolutionStepValue(GRAIN_DIAMETER) = grain_diameter;
199 
201  rGeom[i].FastGetSolutionStepValue(GRAIN_DENSITY) = grain_density;
202 
203  if (mrModelPart.GetNodalSolutionStepVariablesList().Has(REGULARIZATION_COEFFICIENT))
204  rGeom[i].FastGetSolutionStepValue(REGULARIZATION_COEFFICIENT) = regularization_coefficient;
205  }
206  }
207 
208  }
209 
210  KRATOS_CATCH(" ")
211 }; // namespace Kratos
212 
216 
220 
224 
228 
230 std::string Info() const override
231 {
232  return "SetMaterialPropertiesToFluidNodesProcess";
233 }
234 
236 void PrintInfo(std::ostream &rOStream) const override
237 {
238  rOStream << "SetMaterialPropertiesToFluidNodesProcess";
239 }
240 
241 protected:
244 
248 
253 
257 
261 
263 
264 private:
267 
271 
275 
279 
283 
287 
291 
294 
296 // SetMaterialPropertiesToFluidNodesProcess(SetMaterialPropertiesToFluidNodesProcess const& rOther);
297 
299 }
300 ; // Class SetMaterialPropertiesToFluidNodesProcess
301 
303 
306 
310 
312 inline std::istream &operator>>(std::istream &rIStream,
314 
316 inline std::ostream &operator<<(std::ostream &rOStream,
318 {
319  rThis.PrintInfo(rOStream);
320  rOStream << std::endl;
321  rThis.PrintData(rOStream);
322 
323  return rOStream;
324 }
326 
327 } // namespace Kratos.
328 
329 #endif // KRATOS_SET_MATERIAL_PROPERTIES_TO_FLUID_NODES_PROCESS_H_INCLUDED defined
Geometry base class.
Definition: geometry.h:71
SizeType PointsNumber() const
Definition: geometry.h:528
PointerVector< TPointType > PointsArrayType
Definition: geometry.h:118
IndexType Id() const
Definition: indexed_object.h:107
The base class for processes passed to the solution scheme.
Definition: mesher_process.hpp:37
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: mesher_process.hpp:157
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
MeshType::ElementsContainerType ElementsContainerType
Element container. A vector set of Elements with their Id's as key.
Definition: model_part.h:168
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
MeshType::NodesContainerType NodesContainerType
Nodes container. Which is a vector set of nodes with their Id's as key.
Definition: model_part.h:128
MeshType::ElementIterator ElementIterator
Definition: model_part.h:174
VariablesList & GetNodalSolutionStepVariablesList()
Definition: model_part.h:549
static void PartitionedIterators(TVector &rVector, typename TVector::iterator &rBegin, typename TVector::iterator &rEnd)
Generate a partition for an std::vector-like array, providing iterators to the begin and end position...
Definition: openmp_utils.h:179
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
bool Has(TVariableType const &rThisVariable) const
Definition: properties.h:578
Short class definition.
Definition: set_material_properties_to_fluid_nodes_process.hpp:67
virtual ~SetMaterialPropertiesToFluidNodesProcess()
Destructor.
Definition: set_material_properties_to_fluid_nodes_process.hpp:86
std::string Info() const override
Turn back information as a string.
Definition: set_material_properties_to_fluid_nodes_process.hpp:230
void Execute() override
Execute method is used to execute the MesherProcess algorithms.
Definition: set_material_properties_to_fluid_nodes_process.hpp:99
ModelPart & mrModelPart
Definition: set_material_properties_to_fluid_nodes_process.hpp:252
SetMaterialPropertiesToFluidNodesProcess(ModelPart &rModelPart)
Default constructor.
Definition: set_material_properties_to_fluid_nodes_process.hpp:80
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: set_material_properties_to_fluid_nodes_process.hpp:236
void operator()()
Definition: set_material_properties_to_fluid_nodes_process.hpp:90
KRATOS_CLASS_POINTER_DEFINITION(SetMaterialPropertiesToFluidNodesProcess)
Pointer definition of SetMaterialPropertiesToFluidNodesProcess.
bool Has(const VariableData &rThisVariable) const
Definition: variables_list.h:372
#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
ModelPart::NodesContainerType NodesContainerType
Definition: find_conditions_neighbours_process.h:44
std::size_t SizeType
The definition of the size type.
Definition: mortar_classes.h:43
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
ModelPart::ElementsContainerType ElementsContainerType
Definition: clear_contact_conditions_mesher_process.hpp:43
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
GeometryType::PointsArrayType PointsArrayType
Definition: settle_model_structure_process.hpp:48
float viscosity
Definition: edgebased_var.py:8
float density
Definition: face_heat.py:56
bulk_modulus
Definition: script.py:97
integer i
Definition: TensorModule.f:17