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.
calculate_distance_to_skin_process.h
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: Pooyan Dadvand, Ruben Zorrilla
11 //
12 
13 #if !defined(KRATOS_CALCULATE_DISTANCE_TO_SKIN_PROCESS_H_INCLUDED )
14 #define KRATOS_CALCULATE_DISTANCE_TO_SKIN_PROCESS_H_INCLUDED
15 
16 // System includes
17 #include <string>
18 #include <iostream>
19 
20 // External includes
21 
22 // Project includes
26 
27 namespace Kratos
28 {
31 
34 
36 
38 template<std::size_t TDim = 3>
40 {
41 
42 public:
45 
48 
49  // //TODO: These using statements have been included to make the old functions able to compile. It is still pending to update them.
50  // using ConfigurationType = Internals::DistanceSpatialContainersConfigure;
51  // using CellType = OctreeBinaryCell<ConfigurationType>;
52  // using OctreeType = OctreeBinary<CellType>;
53  // using CellNodeDataType = ConfigurationType::cell_node_data_type;
54 
56 
61 
62 
66 
76  ModelPart& rVolumePart,
77  ModelPart& rSkinPart);
78 
90  ModelPart& rVolumePart,
91  ModelPart& rSkinPart,
92  const double RayCastingRelativeTolerance);
93 
104  ModelPart& rVolumePart,
105  ModelPart& rSkinPart,
106  Parameters rParameters);
107 
109  ~CalculateDistanceToSkinProcess() override;
110 
114 
117 
120 
123 
127 
133  void Initialize() override;
134 
143  void CalculateDistances(std::vector<PointerVector<GeometricalObject>>& rIntersectedObjects) override;
144 
153  void CalculateElementalDistances(std::vector<PointerVector<GeometricalObject>> &rIntersectedObjects);
154 
165  double CalculateDistanceToNode(
166  Node &rNode,
167  PointerVector<GeometricalObject> &rIntersectedObjects,
168  const double Epsilon);
169 
174  virtual void InitializeNodalDistances();
175 
183  virtual void CalculateNodalDistances(NodeScalarGetFunctionType& rGetDistanceFunction);
184 
190  virtual void CalculateRayDistances();
191 
195  const Parameters GetDefaultParameters() const override;
196 
200 
202  std::string Info() const override;
203 
205  void PrintInfo(std::ostream& rOStream) const override;
206 
208  void PrintData(std::ostream& rOStream) const override;
209 
211 private:
214 
215 
219 
220  Parameters mSettings;
221 
222  double mRayCastingRelativeTolerance = 1.0e-8;
223 
224  const Variable<double>* mpDistanceVariable = &DISTANCE;
225 
229 
230 
234 
242  double inline CalculatePointDistance(
243  const Element::GeometryType &rIntObjGeom,
244  const Point &rDistancePoint);
245 
249 
250 
254 
255 
259 
260 
262 }; // Class CalculateDistanceToSkinProcess
263 
267 
268 
272 
274 inline std::istream& operator >> (
275  std::istream& rIStream,
277 
279 inline std::ostream& operator << (
280  std::ostream& rOStream,
282 {
283  rThis.PrintInfo(rOStream);
284  rOStream << std::endl;
285  rThis.PrintData(rOStream);
286 
287  return rOStream;
288 }
291 
292 } // namespace Kratos.
293 
294 #endif // KRATOS_CALCULATE_DISTANCE_TO_SKIN_PROCESS_H_INCLUDED defined
std::string Info() const override
Turn back information as a string.
Definition: periodic_interface_process.hpp:93
Applies ray casting to distinguish the color (like in/out) of each node in modelpart.
Definition: apply_ray_casting_process.h:38
This only calculates the distance. Calculating the inside outside should be done by a derived class o...
Definition: calculate_discontinuous_distance_to_skin_process.h:55
Calculates the nodal distances using elemental discontinuous distances.
Definition: calculate_distance_to_skin_process.h:40
CalculateDistanceToSkinProcess & operator=(CalculateDistanceToSkinProcess const &rOther)=delete
Assignment operator.
CalculateDistanceToSkinProcess(CalculateDistanceToSkinProcess const &rOther)=delete
Copy constructor.
typename ApplyRayCastingProcess< TDim >::IntersectionGeometryType IntersectionGeometryType
Types from the ApplyRayCastingProcess.
Definition: calculate_distance_to_skin_process.h:58
CalculateDistanceToSkinProcess()=delete
Default constructor.
typename ApplyRayCastingProcess< TDim >::NodeScalarGetFunctionType NodeScalarGetFunctionType
Definition: calculate_distance_to_skin_process.h:60
typename ApplyRayCastingProcess< TDim >::IntersectionsContainerType IntersectionsContainerType
Definition: calculate_distance_to_skin_process.h:59
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: calculate_distance_to_skin_process.cpp:294
KRATOS_CLASS_POINTER_DEFINITION(CalculateDistanceToSkinProcess)
Pointer definition of CalculateDistanceToSkinProcess.
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: calculate_distance_to_skin_process.cpp:301
Geometry base class.
Definition: geometry.h:71
This class aims to manage meshes for multi-physics simulations.
Definition: model_part.h:77
Node NodeType
Definition: model_part.h:117
This class defines the node.
Definition: node.h:65
This class provides to Kratos a data structure for I/O based on the standard of JSON.
Definition: kratos_parameters.h:59
Point class.
Definition: point.h:59
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
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