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.
geometrical_object_configure.h
Go to the documentation of this file.
1 //
2 // Author: Miquel Santasusana msantasusana@cimne.upc.edu
3 //
4 
5 
6 
7 #if !defined(KRATOS_GEOMETRICAL_OBJECT_CONFIGURE_INCLUDED)
8 #define KRATOS_GEOMETRICAL_OBJECT_CONFIGURE_INCLUDED
9 
10 // System includes
11 #include <string>
12 #include <iostream>
13 #include <cmath>
14 
15 // Kratos includes
17 
18 namespace Kratos
19 {
20 
23 
27 
31 
35 
39 
40 
41 template <std::size_t TDimension>
43 {
44 
45 public:
46 
47  enum {
48  Dimension = TDimension,
49  DIMENSION = TDimension,
50  MAX_LEVEL = 16,
51  MIN_LEVEL = 2
52  };
53 
56 
57 // typedef SpatialSearch SearchType;
58 
59 // typedef SearchType::PointType PointType;
60 // typedef SearchType::ElementsContainerType::ContainerType ContainerType;
61 // typedef SearchType::ElementsContainerType ElementsContainerType;
62 
63 // typedef ContainerType::value_type PointerType;
64 // typedef ContainerType::iterator IteratorType;
65 // typedef ElementsContainerType::iterator ElementIteratorType;
66 //
67 // typedef SearchType::ElementsContainerType::ContainerType ResultContainerType;
68 // // typedef SearchType::ResultDistanceType::ContainerType ResultDistanceType;
69 //
70 // typedef ResultContainerType::iterator ResultIteratorType;
71 // typedef std::vector<double>::iterator DistanceIteratorType;
72 //
73 // typedef ContactPair<PointerType> ContactPairType;
74 // typedef std::vector<ContactPairType> ContainerContactType;
75 // typedef ContainerContactType::iterator IteratorContactType;
76 // typedef ContainerContactType::value_type PointerContactType;
77 
79 
81 
85 
86  typedef ContainerType::value_type PointerType;
87  typedef ContainerType::iterator IteratorType;
88 // typedef ElementsContainerType::iterator ElementIteratorType;
89 
91 // typedef SearchType::ResultDistanceType::ContainerType ResultDistanceType;
92 
93  typedef ResultContainerType::iterator ResultIteratorType;
94  typedef std::vector<double>::iterator DistanceIteratorType;
95 
99 
102  }
103 
107 
108 
112 
113  //******************************************************************************************************************
114 
115  static inline void CalculateBoundingBox(const PointerType& rObject, PointType& rLowPoint, PointType& rHighPoint)
116  {
117  noalias(rHighPoint) = rObject->GetGeometry()[0];
118  noalias(rLowPoint) = rObject->GetGeometry()[0];
119  SphericParticle* p_particle = dynamic_cast<SphericParticle*>(&*rObject);
120  const double& radius = p_particle->GetSearchRadius();
121 
122  for(std::size_t i = 0; i < 3; i++)
123  {
124  rLowPoint[i] += -radius;
125  rHighPoint[i] += radius;
126  }
127  }
128 
129  static inline void CalculateBoundingBox(const PointerType& rObject, PointType& rLowPoint, PointType& rHighPoint, const double& Radius)
130  {
131  noalias(rHighPoint) = rObject->GetGeometry()[0];
132  noalias(rLowPoint) = rObject->GetGeometry()[0];
133 
134  for(std::size_t i = 0; i < 3; i++)
135  {
136  rLowPoint[i] += -Radius;
137  rHighPoint[i] += Radius;
138  }
139  }
140 
141  static inline void CalculateCenter(const PointerType& rObject, PointType& rCenter)
142  {
143  rCenter = rObject->GetGeometry()[0];
144  }
145 
146  //******************************************************************************************************************
147 
148  static inline bool Intersection(const PointerType& rObj_1, const PointerType& rObj_2)
149  {
150  array_1d<double, 3> rObj_2_to_rObj_1;
151  noalias(rObj_2_to_rObj_1) = rObj_1->GetGeometry()[0] - rObj_2->GetGeometry()[0];
152  double distance_2 = inner_prod(rObj_2_to_rObj_1, rObj_2_to_rObj_1);
153 
154  SphericParticle* p_particle1 = dynamic_cast<SphericParticle*>(&*rObj_1);
155  SphericParticle* p_particle2 = dynamic_cast<SphericParticle*>(&*rObj_2);
156  double radius_sum = p_particle1->GetSearchRadius() + p_particle2->GetSearchRadius();
157  bool intersect = floatle((distance_2 - radius_sum * radius_sum),0);
158  return intersect;
159  }
160 
161  static inline bool Intersection(const PointerType& rObj_1, const PointerType& rObj_2, const double& Radius)
162  {
163  array_1d<double, 3> rObj_2_to_rObj_1;
164  noalias(rObj_2_to_rObj_1) = rObj_1->GetGeometry()[0] - rObj_2->GetGeometry()[0];
165  double distance_2 = inner_prod(rObj_2_to_rObj_1, rObj_2_to_rObj_1);
166 
167  SphericParticle* p_particle1 = dynamic_cast<SphericParticle*>(&*rObj_1);
168  SphericParticle* p_particle2 = dynamic_cast<SphericParticle*>(&*rObj_2);
169  double radius_sum = p_particle1->GetSearchRadius() + p_particle2->GetSearchRadius();
170  bool intersect = floatle((distance_2 - radius_sum * radius_sum),0);
171 
172  return intersect;
173  }
174 
175  //******************************************************************************************************************
176 
177  static inline bool IntersectionBox(const PointerType& rObject, const PointType& rLowPoint, const PointType& rHighPoint)
178  {
179  array_1d<double, 3> center_of_particle = rObject->GetGeometry()[0];
180 
181  SphericParticle* p_particle = dynamic_cast<SphericParticle*>(&*rObject);
182  const double& radius = p_particle->GetSearchRadius();
183 
184  bool intersect = (
185  floatle(rLowPoint[0] - radius,center_of_particle[0]) &&
186  floatle(rLowPoint[1] - radius,center_of_particle[1]) &&
187  floatle(rLowPoint[2] - radius,center_of_particle[2]) &&
188  floatge(rHighPoint[0] + radius,center_of_particle[0]) &&
189  floatge(rHighPoint[1] + radius,center_of_particle[1]) &&
190  floatge(rHighPoint[2] + radius,center_of_particle[2]));
191 
192  return intersect;
193  }
194 
195  static inline bool IntersectionBox(const PointerType& rObject, const PointType& rLowPoint, const PointType& rHighPoint, const double& Radius)
196  {
197  array_1d<double, 3> center_of_particle = rObject->GetGeometry()[0];
198 
199  double radius = Radius;//Cambien el radi del objecte de cerca per el gran, aixi no tindria que petar res
200  bool intersect = (
201  floatle(rLowPoint[0] - radius,center_of_particle[0]) &&
202  floatle(rLowPoint[1] - radius,center_of_particle[1]) &&
203  floatle(rLowPoint[2] - radius,center_of_particle[2]) &&
204  floatge(rHighPoint[0] + radius,center_of_particle[0]) &&
205  floatge(rHighPoint[1] + radius,center_of_particle[1]) &&
206  floatge(rHighPoint[2] + radius,center_of_particle[2]));
207  return intersect;
208  }
209 
210  static inline void Distance(const PointerType& rObj_1, const PointerType& rObj_2, double& distance)
211  {
212  array_1d<double, 3> center_of_particle1 = rObj_1->GetGeometry()[0];
213  array_1d<double, 3> center_of_particle2 = rObj_2->GetGeometry()[0];
214 
215  distance = std::sqrt((center_of_particle1[0] - center_of_particle2[0]) * (center_of_particle1[0] - center_of_particle2[0]) +
216  (center_of_particle1[1] - center_of_particle2[1]) * (center_of_particle1[1] - center_of_particle2[1]) +
217  (center_of_particle1[2] - center_of_particle2[2]) * (center_of_particle1[2] - center_of_particle2[2]) );
218  }
219 
220  //******************************************************************************************************************
221 
225 
226 
230 
231 
235 
237  virtual std::string Info() const {return " Spatial Containers Configure for Particles"; }
238 
240  virtual void PrintInfo(std::ostream& rOStream) const {}
241 
243  virtual void PrintData(std::ostream& rOStream) const {}
244 
248 
249 
251 
252 protected:
255 
259 
263 
267 
271 
275 
279 
281 
282 private:
285 
289 
293 
297 
298  static inline bool floateq(double a, double b) {
299  return std::fabs(a - b) < std::numeric_limits<double>::epsilon();
300  }
301 
302  static inline bool floatle(double a, double b) {
303  return std::fabs(a - b) < std::numeric_limits<double>::epsilon() || a < b;
304  }
305 
306  static inline bool floatge(double a, double b) {
307  return std::fabs(a - b) < std::numeric_limits<double>::epsilon() || a > b;
308  }
309 
313 
317 
321 
323  GeometricalConfigure& operator=(GeometricalConfigure const& rOther);
324 
327 
329 
330  }; // Class ParticleConfigure
331 
333 
336 
340 
342  template <std::size_t TDimension>
343  inline std::istream& operator >> (std::istream& rIStream, GeometricalConfigure<TDimension> & rThis){
344  return rIStream;
345  }
346 
348  template <std::size_t TDimension>
349  inline std::ostream& operator << (std::ostream& rOStream, const GeometricalConfigure<TDimension>& rThis){
350  rThis.PrintInfo(rOStream);
351  rOStream << std::endl;
352  rThis.PrintData(rOStream);
353 
354  return rOStream;
355  }
356 
358 
359 } // namespace Kratos.
360 #endif /* KRATOS_GEOMETRICAL_OBJECT_CONFIGURE_H */
Definition: geometrical_object_configure.h:43
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: geometrical_object_configure.h:243
KRATOS_CLASS_POINTER_DEFINITION(GeometricalConfigure)
Pointer definition of SpatialContainersConfigure.
virtual ~GeometricalConfigure()
Definition: geometrical_object_configure.h:101
static void CalculateBoundingBox(const PointerType &rObject, PointType &rLowPoint, PointType &rHighPoint)
Definition: geometrical_object_configure.h:115
static void CalculateCenter(const PointerType &rObject, PointType &rCenter)
Definition: geometrical_object_configure.h:141
static bool Intersection(const PointerType &rObj_1, const PointerType &rObj_2)
Definition: geometrical_object_configure.h:148
ResultContainerType::iterator ResultIteratorType
Definition: geometrical_object_configure.h:93
static void Distance(const PointerType &rObj_1, const PointerType &rObj_2, double &distance)
Definition: geometrical_object_configure.h:210
virtual std::string Info() const
Turn back information as a string.
Definition: geometrical_object_configure.h:237
PointerVectorSet< GeometricalObject, IndexedObject >::ContainerType ResultContainerType
Definition: geometrical_object_configure.h:90
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: geometrical_object_configure.h:240
ContainerType::iterator IteratorType
Definition: geometrical_object_configure.h:87
GeometricalConfigure()
Definition: geometrical_object_configure.h:100
static bool IntersectionBox(const PointerType &rObject, const PointType &rLowPoint, const PointType &rHighPoint)
Definition: geometrical_object_configure.h:177
SpatialSearch SearchType
Definition: geometrical_object_configure.h:80
PointerVectorSet< GeometricalObject, IndexedObject > ElementsContainerType
Definition: geometrical_object_configure.h:84
ContainerType::value_type PointerType
Definition: geometrical_object_configure.h:86
static bool Intersection(const PointerType &rObj_1, const PointerType &rObj_2, const double &Radius)
Definition: geometrical_object_configure.h:161
SearchType::PointType PointType
Definition: geometrical_object_configure.h:82
static void CalculateBoundingBox(const PointerType &rObject, PointType &rLowPoint, PointType &rHighPoint, const double &Radius)
Definition: geometrical_object_configure.h:129
PointerVectorSet< GeometricalObject, IndexedObject >::ContainerType ContainerType
Definition: geometrical_object_configure.h:83
@ Dimension
Definition: geometrical_object_configure.h:48
@ MIN_LEVEL
Definition: geometrical_object_configure.h:51
@ DIMENSION
Definition: geometrical_object_configure.h:49
@ MAX_LEVEL
Definition: geometrical_object_configure.h:50
std::vector< double >::iterator DistanceIteratorType
Definition: geometrical_object_configure.h:94
static bool IntersectionBox(const PointerType &rObject, const PointType &rLowPoint, const PointType &rHighPoint, const double &Radius)
Definition: geometrical_object_configure.h:195
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
Definition: spheric_particle.h:31
virtual double GetSearchRadius()
Definition: spheric_particle.cpp:2201
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
TExpression1Type::data_type inner_prod(AMatrix::MatrixExpression< TExpression1Type, TCategory1 > const &First, AMatrix::MatrixExpression< TExpression2Type, TCategory2 > const &Second)
Definition: amatrix_interface.h:592
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
T & noalias(T &TheMatrix)
Definition: amatrix_interface.h:484
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
a
Definition: generate_stokes_twofluid_element.py:77
b
Definition: generate_total_lagrangian_mixed_volumetric_strain_element.py:31
float radius
Definition: mesh_to_mdpa_converter.py:18
distance_2
Definition: sp_statistics.py:71
integer i
Definition: TensorModule.f:17