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.
find_nodal_h_for_rigid_walls_process.hpp
Go to the documentation of this file.
1 //
2 // Project Name: KratosPfemFluidDynamicsApplication $
3 // Created by: $Author: AFranci $
4 // Last modified by: $Co-Author: $
5 // Date: $Date: October 2019 $
6 // Revision: $Revision: 0.0 $
7 //
8 //
9 
10 #if !defined(KRATOS_FIND_NODAL_H_FOR_RIGID_WALLS_PROCESS_H_INCLUDED)
11 #define KRATOS_FIND_NODAL_H_FOR_RIGID_WALLS_PROCESS_H_INCLUDED
12 
13 // External includes
14 
15 // System includes
16 
17 // Project includes
20 
21 #include "includes/model_part.h"
24 
27 
29 // Data:
30 // Flags: (checked)
31 // (set)
32 // (modified)
33 // (reset)
34 //(set):=(set in this process)
35 
36 namespace Kratos
37 {
38 
41 
43 
48  : public Process
49  {
50  public:
53 
56 
61  typedef std::size_t SizeType;
65 
68  : mrModelPart(rModelPart)
69  {
70  KRATOS_INFO("FindNodalHForRigidWallsProcess") << " activated " << std::endl;
71  }
72 
75 
79 
81  void operator()()
82  {
83  Execute();
84  }
85 
89 
90  void Execute() override
91  {
93  // Check if variables are available
94 
95  // initialize to zero
96  for (ModelPart::NodesContainerType::iterator i_node = mrModelPart.NodesBegin(); i_node != mrModelPart.NodesEnd(); i_node++)
97  {
98  if (i_node->Is(RIGID))
99  {
100  i_node->FastGetSolutionStepValue(NODAL_H_WALL) = 0;
101  }
102  }
103 
104  for (IndexType i = 0; i < mrModelPart.Elements().size(); ++i)
105  {
106  auto it_element = mrModelPart.ElementsBegin() + i;
107  auto &r_geom = it_element->GetGeometry();
108  const SizeType number_of_nodes = r_geom.size();
109 
110  for (IndexType k = 0; k < number_of_nodes - 1; ++k)
111  {
112  if (r_geom[k].Is(RIGID))
113  {
114  double &r_h1 = r_geom[k].FastGetSolutionStepValue(NODAL_H_WALL);
115  for (IndexType l = k + 1; l < number_of_nodes; ++l)
116  {
117  if (r_geom[l].Is(RIGID))
118  {
119  double hedge = norm_2(r_geom[l].Coordinates() - r_geom[k].Coordinates());
120  double &r_h2 = r_geom[l].FastGetSolutionStepValue(NODAL_H_WALL);
121  // Get minimum between the existent value and the considered edge length
122  r_geom[k].FastGetSolutionStepValue(NODAL_H_WALL) = std::max(r_h1, hedge);
123  r_geom[l].FastGetSolutionStepValue(NODAL_H_WALL) = std::max(r_h2, hedge);
124  }
125  }
126  }
127  }
128  }
129  // Synchronize between processes
130  // mrModelPart.GetCommunicator().SynchronizeCurrentDataToMin(NODAL_H_WALL);
131  KRATOS_CATCH("")
132  }
133 
137 
141 
145 
147  std::string Info() const override
148  {
149  return "FindNodalHForRigidWallsProcess";
150  }
151 
153  void PrintInfo(std::ostream &rOStream) const override
154  {
155  rOStream << "FindNodalHForRigidWallsProcess";
156  }
157 
161 
163 
164  private:
167 
171  ModelPart &mrModelPart;
172 
176 
180 
184 
188 
192 
195 
197 
199  // Process(Process const& rOther);
200 
202 
203  }; // Class Process
204 
206 
209 
213 
215  inline std::istream &operator>>(std::istream &rIStream,
217 
219  inline std::ostream &operator<<(std::ostream &rOStream,
220  const FindNodalHForRigidWallsProcess &rThis)
221  {
222  rThis.PrintInfo(rOStream);
223  rOStream << std::endl;
224  rThis.PrintData(rOStream);
225 
226  return rOStream;
227  }
229 
230 } // namespace Kratos.
231 
232 #endif // KRATOS_SET_INLET_PROCESS_H_INCLUDED defined
Base class for all Conditions.
Definition: condition.h:59
Refine Mesh Elements Process 2D and 3D.
Definition: find_nodal_h_for_rigid_walls_process.hpp:49
ModelPart::PropertiesType PropertiesType
Definition: find_nodal_h_for_rigid_walls_process.hpp:59
ModelPart::ConditionType ConditionType
Definition: find_nodal_h_for_rigid_walls_process.hpp:58
std::string Info() const override
Turn back information as a string.
Definition: find_nodal_h_for_rigid_walls_process.hpp:147
ConditionType::GeometryType GeometryType
Definition: find_nodal_h_for_rigid_walls_process.hpp:60
ModelPart::NodeType NodeType
Definition: find_nodal_h_for_rigid_walls_process.hpp:57
virtual ~FindNodalHForRigidWallsProcess()
Destructor.
Definition: find_nodal_h_for_rigid_walls_process.hpp:74
FindNodalHForRigidWallsProcess(ModelPart &rModelPart)
Default constructor.
Definition: find_nodal_h_for_rigid_walls_process.hpp:67
KRATOS_CLASS_POINTER_DEFINITION(FindNodalHForRigidWallsProcess)
Pointer definition of Process.
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: find_nodal_h_for_rigid_walls_process.hpp:153
void Execute() override
Execute method is used to execute the Process algorithms.
Definition: find_nodal_h_for_rigid_walls_process.hpp:90
std::size_t SizeType
Definition: find_nodal_h_for_rigid_walls_process.hpp:61
void operator()()
This operator is provided to call the process as a function and simply calls the Execute method.
Definition: find_nodal_h_for_rigid_walls_process.hpp:81
std::size_t IndexType
Definition: flags.h:74
bool Is(Flags const &rOther) const
Definition: flags.h:274
Geometry base class.
Definition: geometry.h:71
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
ElementIterator ElementsBegin(IndexType ThisIndex=0)
Definition: model_part.h:1169
NodeIterator NodesBegin(IndexType ThisIndex=0)
Definition: model_part.h:487
ElementsContainerType & Elements(IndexType ThisIndex=0)
Definition: model_part.h:1189
NodeIterator NodesEnd(IndexType ThisIndex=0)
Definition: model_part.h:497
This class defines the node.
Definition: node.h:65
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
Properties encapsulates data shared by different Elements or Conditions. It can store any type of dat...
Definition: properties.h:69
#define KRATOS_CATCH(MoreInfo)
Definition: define.h:110
#define KRATOS_TRY
Definition: define.h:109
#define KRATOS_INFO(label)
Definition: logger.h:250
static double max(double a, double b)
Definition: GeometryFunctions.h:79
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
TExpressionType::data_type norm_2(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:625
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 k
Definition: quadrature.py:595
integer i
Definition: TensorModule.f:17
integer l
Definition: TensorModule.f:17