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_from_fluid_to_rigid_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_FROM_FLUID_TO_RIGID_NODES_PROCESS_H_INCLUDED)
11 #define KRATOS_SET_MATERIAL_PROPERTIES_FROM_FLUID_TO_RIGID_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  : mrRigidModelPart(rRigidModelPart),
82  mrFluidModelPart(rFluidModelPart)
83  {
84  }
85 
88  {
89  }
90 
91  void operator()()
92  {
93  Execute();
94  }
95 
99 
100  void Execute() override
101  {
102  KRATOS_TRY
103 
104  double density = 0;
105  double bulk_modulus = 0;
106  double viscosity = 0;
107  double flow_index = 1;
108  double yield_shear = 0;
109  double adaptive_exponent = 0;
110  double static_friction = 0;
111  double dynamic_friction = 0;
112  double inertial_number_zero = 0;
113  double grain_diameter = 0;
114  double grain_density = 0;
115  double regularization_coefficient = 0;
116  double friction_angle = 0;
117  double cohesion = 0;
118 
119 #pragma omp parallel
120  {
121 
122  ModelPart::ElementIterator ElemBegin;
125  for (ModelPart::ElementIterator itElem = ElemBegin; itElem != ElemEnd; ++itElem)
126  {
127  ModelPart::PropertiesType &elemProperties = itElem->GetProperties();
128 
129  density = elemProperties[DENSITY];
130  bulk_modulus = elemProperties[BULK_MODULUS];
131  viscosity = elemProperties[DYNAMIC_VISCOSITY];
132 
133  if (elemProperties.Has(YIELD_SHEAR)) // Bingham model
134  {
135  flow_index = elemProperties[FLOW_INDEX];
136  yield_shear = elemProperties[YIELD_SHEAR];
137  adaptive_exponent = elemProperties[ADAPTIVE_EXPONENT];
138  }
139  else if (elemProperties.Has(INTERNAL_FRICTION_ANGLE)) // Frictional Viscoplastic model
140  {
141  friction_angle = elemProperties[INTERNAL_FRICTION_ANGLE];
142  cohesion = elemProperties[COHESION];
143  adaptive_exponent = elemProperties[ADAPTIVE_EXPONENT];
144  }
145  else if (elemProperties.Has(STATIC_FRICTION)) // Mu(I)-rheology
146  {
147  static_friction = elemProperties[STATIC_FRICTION];
148  dynamic_friction = elemProperties[DYNAMIC_FRICTION];
149  inertial_number_zero = elemProperties[INERTIAL_NUMBER_ZERO];
150  grain_diameter = elemProperties[GRAIN_DIAMETER];
151  grain_density = elemProperties[GRAIN_DENSITY];
152  regularization_coefficient = elemProperties[REGULARIZATION_COEFFICIENT];
153  }
154  break;
155  }
156  }
157 #pragma omp parallel
158  {
159 
160  ModelPart::NodeIterator NodeBegin;
161  ModelPart::NodeIterator NodeEnd;
163  for (ModelPart::NodeIterator iNode = NodeBegin; iNode != NodeEnd; ++iNode)
164  {
165 
167  iNode->FastGetSolutionStepValue(BULK_MODULUS) = bulk_modulus;
168 
170  iNode->FastGetSolutionStepValue(DENSITY) = density;
171 
172  if (mrFluidModelPart.GetNodalSolutionStepVariablesList().Has(DYNAMIC_VISCOSITY))
173  iNode->FastGetSolutionStepValue(DYNAMIC_VISCOSITY) = viscosity;
174 
176  iNode->FastGetSolutionStepValue(YIELD_SHEAR) = yield_shear;
177 
179  iNode->FastGetSolutionStepValue(FLOW_INDEX) = flow_index;
180 
181  if (mrFluidModelPart.GetNodalSolutionStepVariablesList().Has(ADAPTIVE_EXPONENT))
182  iNode->FastGetSolutionStepValue(ADAPTIVE_EXPONENT) = adaptive_exponent;
183 
184  if (mrFluidModelPart.GetNodalSolutionStepVariablesList().Has(INTERNAL_FRICTION_ANGLE))
185  iNode->FastGetSolutionStepValue(INTERNAL_FRICTION_ANGLE) = friction_angle;
186 
188  iNode->FastGetSolutionStepValue(COHESION) = cohesion;
189 
190  if (mrFluidModelPart.GetNodalSolutionStepVariablesList().Has(ADAPTIVE_EXPONENT))
191  iNode->FastGetSolutionStepValue(ADAPTIVE_EXPONENT) = adaptive_exponent;
192 
194  iNode->FastGetSolutionStepValue(STATIC_FRICTION) = static_friction;
195 
197  iNode->FastGetSolutionStepValue(DYNAMIC_FRICTION) = dynamic_friction;
198 
199  if (mrFluidModelPart.GetNodalSolutionStepVariablesList().Has(INERTIAL_NUMBER_ZERO))
200  iNode->FastGetSolutionStepValue(INERTIAL_NUMBER_ZERO) = inertial_number_zero;
201 
203  iNode->FastGetSolutionStepValue(GRAIN_DIAMETER) = grain_diameter;
204 
206  iNode->FastGetSolutionStepValue(GRAIN_DENSITY) = grain_density;
207 
208  if (mrFluidModelPart.GetNodalSolutionStepVariablesList().Has(REGULARIZATION_COEFFICIENT))
209  iNode->FastGetSolutionStepValue(REGULARIZATION_COEFFICIENT) = regularization_coefficient;
210  }
211  }
212 
213  KRATOS_CATCH(" ")
214  };
215 
219 
223 
227 
231 
233  std::string Info() const override
234  {
235  return "SetMaterialPropertiesFromFluidToRigidNodesProcess";
236  }
237 
239  void PrintInfo(std::ostream &rOStream) const override
240  {
241  rOStream << "SetMaterialPropertiesFromFluidToRigidNodesProcess";
242  }
243 
244  protected:
247 
251 
257 
261 
265 
267 
268  private:
271 
275 
279 
283 
287 
291 
295 
298 
300  // SetMaterialPropertiesFromFluidToRigidNodesProcess(SetMaterialPropertiesFromFluidToRigidNodesProcess const& rOther);
301 
303 
304  }; // namespace Kratos
305 
307 
310 
314 
316  inline std::istream &operator>>(std::istream &rIStream,
318 
320  inline std::ostream &operator<<(std::ostream &rOStream,
322  {
323  rThis.PrintInfo(rOStream);
324  rOStream << std::endl;
325  rThis.PrintData(rOStream);
326 
327  return rOStream;
328  }
330 
331 } // namespace Kratos.
332 
333 #endif // KRATOS_SET_MATERIAL_PROPERTIES_FROM_FLUID_TO_RIGID_NODES_PROCESS_H_INCLUDED defined
PointerVector< TPointType > PointsArrayType
Definition: geometry.h:118
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
MeshType::NodeIterator NodeIterator
Definition: model_part.h:134
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
NodesContainerType & Nodes(IndexType ThisIndex=0)
Definition: model_part.h:507
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_from_fluid_to_rigid_nodes_process.hpp:67
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: set_material_properties_from_fluid_to_rigid_nodes_process.hpp:239
SetMaterialPropertiesFromFluidToRigidNodesProcess(ModelPart &rRigidModelPart, ModelPart &rFluidModelPart)
Default constructor.
Definition: set_material_properties_from_fluid_to_rigid_nodes_process.hpp:80
std::string Info() const override
Turn back information as a string.
Definition: set_material_properties_from_fluid_to_rigid_nodes_process.hpp:233
void operator()()
Definition: set_material_properties_from_fluid_to_rigid_nodes_process.hpp:91
ModelPart & mrRigidModelPart
Definition: set_material_properties_from_fluid_to_rigid_nodes_process.hpp:255
KRATOS_CLASS_POINTER_DEFINITION(SetMaterialPropertiesFromFluidToRigidNodesProcess)
Pointer definition of SetMaterialPropertiesFromFluidToRigidNodesProcess.
ModelPart & mrFluidModelPart
Definition: set_material_properties_from_fluid_to_rigid_nodes_process.hpp:256
virtual ~SetMaterialPropertiesFromFluidToRigidNodesProcess()
Destructor.
Definition: set_material_properties_from_fluid_to_rigid_nodes_process.hpp:87
void Execute() override
Execute method is used to execute the MesherProcess algorithms.
Definition: set_material_properties_from_fluid_to_rigid_nodes_process.hpp:100
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