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.
point_configure.h
Go to the documentation of this file.
1 // Project Name: Kratos
2 // Last Modified by: $Author: Nelson Lafontaine $
3 // Date: $Date: 2012-05-24 $
4 // Revision: $Revision: 1.0 $
5 //
6 
7 
8 #if !defined(KRATOS_POINT_CONFIGURE)
9 #define KRATOS_POINT_CONFIGURE
10 
11 
12 
13 // System includes
14 #include <string>
15 #include <iostream>
16 #include <cmath>
18 
19 // Kratos includes
20 #include "includes/variables.h"
22 //#include "GeometryFunctions.h"
23 
24 namespace Kratos
25 {
26 
29 
33 
37 
41 
45 
46 
47 template <std::size_t TDimension>
49 public:
50 
51 enum {Dimension = TDimension,
52  DIMENSION = TDimension,
53  MAX_LEVEL = 16,
54  MIN_LEVEL = 2
55 };
56 
58 
62 
63 typedef ContainerType::value_type PointerType;
64 typedef ContainerType::iterator IteratorType;
65 
67 
68 
69 typedef ResultContainerType::iterator ResultIteratorType;
70 typedef std::vector<double>::iterator DistanceIteratorType;
71 
74 
78 
80 virtual ~PointConfigure(){}
81 
85 
86 
90 
91 
92 //******************************************************************************************************************
93 //******************************************************************************************************************
94 
95 static inline void CalculateBoundingBox(const PointerType& r_p_point, PointType& r_low_point, PointType& r_high_point){
96  r_high_point = r_low_point = *r_p_point;
97 }
98 
99 //******************************************************************************************************************
100 //******************************************************************************************************************
101 
102 static inline void CalculateBoundingBox(const PointerType& r_p_point, PointType& r_low_point, PointType& r_high_point, const double& radius){
103  r_high_point = r_low_point = *r_p_point;
104 
105  for (std::size_t i = 0; i < 3; ++i){
106  r_low_point[i] -= radius;
107  r_high_point[i] += radius;
108  }
109 }
110 
111 //******************************************************************************************************************
112 //******************************************************************************************************************
113 
114 static inline void CalculateCenter(const PointerType& r_p_point, PointType& rCenter){
115  rCenter = *r_p_point;
116 }
117 
118 //******************************************************************************************************************
119 //******************************************************************************************************************
120 
121 static inline bool Intersection(const PointerType& r_p_point_1, const PointerType& r_p_point_2){
122  return false;
123 }
124 
125 //******************************************************************************************************************
126 //******************************************************************************************************************
127 
128 static inline bool Intersection(const PointerType& r_p_point_1, const PointerType& r_p_point_2, const double& radius){
129  array_1d<double, 3> node_2_to_1 = *r_p_point_1 - *r_p_point_2;
130  double distance;
131  Distance(r_p_point_1, r_p_point_2, distance);
132 
133  bool intersect = (distance - radius) <= 0;
134 
135  return intersect;
136 }
137 
138 //******************************************************************************************************************
139 //******************************************************************************************************************
140 
141 static inline bool IntersectionBox(const PointerType& r_p_point, const PointType& r_low_point, const PointType& r_high_point)
142 {
143  array_1d<double, 3> center = *r_p_point;
144 
145  bool intersect = (r_low_point[0] <= center[0] && r_low_point[1] <= center[1] && r_low_point[2] <= center[2] &&
146  r_high_point[0] >= center[0] && r_high_point[1] >= center[1] && r_high_point[2] >= center[2]);
147 
148  return intersect;
149 
150 }
151 
152 //******************************************************************************************************************
153 //******************************************************************************************************************
154 
155 static inline bool IntersectionBox(const PointerType& r_p_point, const PointType& r_low_point, const PointType& r_high_point, const double& radius)
156 {
157  array_1d<double, 3> center = *r_p_point;
158 
159  bool intersect = (r_low_point[0] - radius <= center[0] && r_low_point[1] - radius <= center[1] && r_low_point[2] - radius <= center[2] &&
160  r_high_point[0] + radius >= center[0] && r_high_point[1] + radius >= center[1] && r_high_point[2] + radius >= center[2]);
161 
162  return intersect;
163 }
164 
165 //******************************************************************************************************************
166 //******************************************************************************************************************
167 
168 static inline void Distance(const PointerType& r_p_point_1, const PointerType& r_p_point_2, double& distance)
169 {
170  array_1d<double, 3> center_1 = *r_p_point_1;
171  array_1d<double, 3> center_2 = *r_p_point_2;
172 
173  distance = sqrt((center_1[0] - center_2[0]) * (center_1[0] - center_2[0]) +
174  (center_1[1] - center_2[1]) * (center_1[1] - center_2[1]) +
175  (center_1[2] - center_2[2]) * (center_1[2] - center_2[2]));
176 }
177 
178 //******************************************************************************************************************
179 //******************************************************************************************************************
180 
184 
185 
189 
190 
194 
196 virtual std::string Info() const {return " Spatial Containers Configure for Particles"; }
197 
199 virtual void PrintInfo(std::ostream& rOStream) const {}
200 
202 virtual void PrintData(std::ostream& rOStream) const {}
203 
207 
208 
210 
211 protected:
214 
218 
222 
226 
230 
234 
238 
240 
241 private:
244 
248 
252 
256 
260 
264 
268 
270 PointConfigure& operator=(PointConfigure const& rOther);
271 
273 PointConfigure(PointConfigure const& rOther);
274 
276 
277 }; // Class ParticleConfigure
278 
280 
283 
287 
289  template <std::size_t TDimension>
290  inline std::istream& operator >> (std::istream& rIStream, PointConfigure<TDimension> & rThis){
291  return rIStream;
292  }
293 
295  template <std::size_t TDimension>
296  inline std::ostream& operator << (std::ostream& rOStream, const PointConfigure<TDimension>& rThis){
297  rThis.PrintInfo(rOStream);
298  rOStream << std::endl;
299  rThis.PrintData(rOStream);
300 
301  return rOStream;
302  }
303 
305 
306 } // namespace Kratos.
307 #endif /* KRATOS_POINT_CONFIGURE */
Definition: point_configure.h:48
static void Distance(const PointerType &r_p_point_1, const PointerType &r_p_point_2, double &distance)
Definition: point_configure.h:168
static bool IntersectionBox(const PointerType &r_p_point, const PointType &r_low_point, const PointType &r_high_point)
Definition: point_configure.h:141
ResultContainerType::iterator ResultIteratorType
Definition: point_configure.h:69
static void CalculateBoundingBox(const PointerType &r_p_point, PointType &r_low_point, PointType &r_high_point, const double &radius)
Definition: point_configure.h:102
PointConfigure()
Definition: point_configure.h:79
ContainerType::iterator IteratorType
Definition: point_configure.h:64
PointerVectorSet< Point, IndexedObject > PointsContainerType
Definition: point_configure.h:61
PointerVectorSet< Point, IndexedObject >::ContainerType ResultContainerType
Definition: point_configure.h:66
static void CalculateBoundingBox(const PointerType &r_p_point, PointType &r_low_point, PointType &r_high_point)
Definition: point_configure.h:95
virtual std::string Info() const
Turn back information as a string.
Definition: point_configure.h:196
SearchType::PointType PointType
Definition: point_configure.h:59
SpatialSearch SearchType
Definition: point_configure.h:57
static bool IntersectionBox(const PointerType &r_p_point, const PointType &r_low_point, const PointType &r_high_point, const double &radius)
Definition: point_configure.h:155
KRATOS_CLASS_POINTER_DEFINITION(PointConfigure)
Pointer definition of SpatialContainersConfigure.
static bool Intersection(const PointerType &r_p_point_1, const PointerType &r_p_point_2)
Definition: point_configure.h:121
static void CalculateCenter(const PointerType &r_p_point, PointType &rCenter)
Definition: point_configure.h:114
PointerVectorSet< Point, IndexedObject >::ContainerType ContainerType
Definition: point_configure.h:60
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: point_configure.h:199
static bool Intersection(const PointerType &r_p_point_1, const PointerType &r_p_point_2, const double &radius)
Definition: point_configure.h:128
@ DIMENSION
Definition: point_configure.h:52
@ MIN_LEVEL
Definition: point_configure.h:54
@ Dimension
Definition: point_configure.h:51
@ MAX_LEVEL
Definition: point_configure.h:53
std::vector< double >::iterator DistanceIteratorType
Definition: point_configure.h:70
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: point_configure.h:202
virtual ~PointConfigure()
Definition: point_configure.h:80
ContainerType::value_type PointerType
Definition: point_configure.h:63
Point class.
Definition: point.h:59
A sorted associative container similar to an STL set, but uses a vector to store pointers to its data...
Definition: pointer_vector_set.h:72
TContainerType ContainerType
Definition: pointer_vector_set.h:90
This class is used to search for elements, conditions and nodes in a given model part.
Definition: spatial_search.h:50
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
float radius
Definition: mesh_to_mdpa_converter.py:18
integer i
Definition: TensorModule.f:17