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_3d.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: Riccardo Rossi
11 // Janosch Stascheit
12 // Felix Nagel
13 // contributors: Hoang Giang Bui
14 // Josep Maria Carbonell
15 //
16 
17 #if !defined(KRATOS_POINT_3D_H_INCLUDED )
18 #define KRATOS_POINT_3D_H_INCLUDED
19 
20 // System includes
21 
22 // External includes
23 
24 // Project includes
25 #include "geometries/geometry.h"
26 
27 namespace Kratos
28 {
29 
32 
36 
40 
44 
48 
51 template<class TPointType>
52 class Point3D : public Geometry<TPointType>
53 {
54 public:
58 
61 
64 
68 
73 
76  typedef TPointType PointType;
77 
82  typedef typename BaseType::IndexType IndexType;
83 
84 
89  typedef typename BaseType::SizeType SizeType;
90 
95 
101 
108 
114 
119 
124 
130 
136 
140 
145 
149 
150 // Point3D(const PointType& FirstPoint)
151 // : BaseType(PointsArrayType(), &msGeometryData)
152 // {
153 // BaseType::Points().push_back(typename PointType::Pointer(new PointType(FirstPoint)));
154 // }
155 
156  Point3D(typename PointType::Pointer pFirstPoint)
157  : BaseType(PointsArrayType(), &msGeometryData)
158  {
159  BaseType::Points().push_back(pFirstPoint);
160  }
161 
162  Point3D(const PointsArrayType& ThisPoints)
163  : BaseType(ThisPoints, &msGeometryData)
164  {
165  if( BaseType::PointsNumber() != 1)
166  KRATOS_ERROR << "Invalid points number. Expected 2, given " << BaseType::PointsNumber() << std::endl;
167  }
168 
170  explicit Point3D(
171  const IndexType GeometryId,
172  const PointsArrayType& rThisPoints
173  ) : BaseType(GeometryId, rThisPoints, &msGeometryData)
174  {
175  KRATOS_ERROR_IF( this->PointsNumber() != 1 ) << "Invalid points number. Expected 1, given " << this->PointsNumber() << std::endl;
176  }
177 
179  explicit Point3D(
180  const std::string& rGeometryName,
181  const PointsArrayType& rThisPoints
182  ) : BaseType(rGeometryName, rThisPoints, &msGeometryData)
183  {
184  KRATOS_ERROR_IF(this->PointsNumber() != 1) << "Invalid points number. Expected 1, given " << this->PointsNumber() << std::endl;
185  }
186 
195  Point3D(Point3D const& rOther)
196  : BaseType(rOther)
197  {
198  }
199 
200 
212  template<class TOtherPointType> Point3D(Point3D<TOtherPointType> const& rOther)
213  : BaseType(rOther)
214  {
215  }
216 
218  ~Point3D() override {}
219 
221  {
223  }
225  {
227  }
228 
232 
243  Point3D& operator=(const Point3D& rOther)
244  {
245  BaseType::operator=(rOther);
246 
247  return *this;
248  }
249 
260  template<class TOtherPointType>
262  {
263  BaseType::operator=(rOther);
264 
265  return *this;
266  }
267 
271 
278  typename BaseType::Pointer Create(
279  const IndexType NewGeometryId,
280  PointsArrayType const& rThisPoints
281  ) const override
282  {
283  return typename BaseType::Pointer( new Point3D( NewGeometryId, rThisPoints ) );
284  }
285 
292  typename BaseType::Pointer Create(
293  const IndexType NewGeometryId,
294  const BaseType& rGeometry
295  ) const override
296  {
297  auto p_geometry = typename BaseType::Pointer( new Point3D( NewGeometryId, rGeometry.Points() ) );
298  p_geometry->SetData(rGeometry.GetData());
299  return p_geometry;
300  }
301 
302 // /**
303 // * @brief Lumping factors for the calculation of the lumped mass matrix
304 // * @param rResult Vector containing the lumping factors
305 // * @param LumpingMethod The lumping method considered. The three methods available are:
306 // * - The row sum method
307 // * - Diagonal scaling
308 // * - Evaluation of M using a quadrature involving only the nodal points and thus automatically yielding a diagonal matrix for standard element shape function
309 // */
310 // Vector& LumpingFactors(
311 // Vector& rResult,
312 // const typename BaseType::LumpingMethods LumpingMethod = BaseType::LumpingMethods::ROW_SUM
313 // ) const override
314 // {
315 // }
316 
320 
333  double Length() const override
334  {
335  return 0.00;
336  }
337 
349  double Area() const override
350  {
351  return 0.00;
352  }
353 
354 
365  double DomainSize() const override
366  {
367  return 0.00;
368  }
369 
373 
388  // virtual JacobiansType& Jacobian(JacobiansType& rResult, IntegrationMethod ThisMethod) const
389  // {
390  // }
391 
409  // virtual Matrix& Jacobian(Matrix& rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const
410  // {
411  // }
412 
424 //
436  // virtual Vector& DeterminantOfJacobian(Vector& rResult, IntegrationMethod ThisMethod) const
437  // {
438  // }
439 
458  // virtual double DeterminantOfJacobian(IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const
459  // {
460  // }
461 
474 // virtual double DeterminantOfJacobian(const CoordinatesArrayType& rPoint) const
475 // {
476 // }
477 
478 
493  // virtual JacobiansType& InverseOfJacobian(JacobiansType& rResult, IntegrationMethod ThisMethod) const
494  // {
495  // }
496 
514  // virtual Matrix& InverseOfJacobian(Matrix& rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const
515  // {
516  // }
517 
529  // virtual Matrix& InverseOfJacobian(Matrix& rResult, const CoordinatesArrayType& rPoint) const
530  // {
531  // }
532 
533 
537  SizeType EdgesNumber() const override
538  {
539  return 1;
540  }
541 
542  SizeType FacesNumber() const override
543  {
544  return 0;
545  }
546 
550 
551 // virtual double ShapeFunctionValue(IndexType ShapeFunctionIndex,
552 // const CoordinatesArrayType& rPoint) const
553 // {
554 // KRATOS_ERROR << "This method is not implemented yet!" << *this << std::endl;
555 // return 0;
556 // }
557 //
558  /*
559 
560 
561  virtual void BoundingBox(Bounding_Box<TPointType>& rResult) const
562  {
563  this->BoundingBox(rResult.mLowPoint, rResult.mHighPoint);
564  }
565 
566 
567  virtual void BoundingBox(TPointType& rLowPoint, TPointType& rHighPoint) const
568  {
569 
570  array_1d<double,3 > point = ZeroVector(3);
571 
572  point[0] = BaseType::GetPoint(0).X();
573  point[1] = BaseType::GetPoint(0).Y();
574  point[2] = BaseType::GetPoint(0).Z();
575 
577  rLowPoint.push_back(point);
578  rHighPoint.push_back(point);
579 
580  }
581 
582  */
583 
587 
594  std::string Info() const override
595  {
596  return "a point in 3D space";
597  }
598 
605  void PrintInfo(std::ostream& rOStream) const override
606  {
607  rOStream << "a point in 3D space";
608  }
609 
618  void PrintData(std::ostream& rOStream) const override
619  {
620  rOStream << "a point in 3D space";
621  }
622 
623 
627 
628 
630 
631 protected:
634 
635 
639 
640 
644 
645 
649 
650 
654 
655 
659 
660 
664 
665 
667 
668 private:
671 
672  static const GeometryData msGeometryData;
673 
674  static const GeometryDimension msGeometryDimension;
675 
679 
680 
684 
685 
689 
690 // static Matrix CalculateShapeFunctionsIntegrationPointsValues(typename BaseType::IntegrationMethod ThisMethod)
691 // {
692 // }
693 
694  // static ShapeFunctionsGradientsType CalculateShapeFunctionsIntegrationPointsLocalGradients(typename BaseType::IntegrationMethod ThisMethod)
695  // {
696  // }
697 
698  static const IntegrationPointsContainerType AllIntegrationPoints()
699  {
700  IntegrationPointsContainerType integration_points;
701  return integration_points;
702  }
703 
704  static const ShapeFunctionsValuesContainerType AllShapeFunctionsValues()
705  {
706  ShapeFunctionsValuesContainerType shape_functions_value;
707  return shape_functions_value;
708  }
709  static const ShapeFunctionsLocalGradientsContainerType AllShapeFunctionsLocalGradients()
710  {
711  ShapeFunctionsLocalGradientsContainerType shape_functions_local_gradients;
712  return shape_functions_local_gradients;
713  }
714 
718 
719 
723 
724 
728 
729  template<class TOtherPointType> friend class Point3D;
730 
734 
735  friend class Serializer;
736 
737  void save(Serializer& rSerializer) const override
738  {
740  }
741 
742  void load(Serializer& rSerializer) override
743  {
745  }
746 
747  // Default constructor needed for serialization only
748  Point3D():BaseType( PointsArrayType(), &msGeometryData ) {}
749 
753 
754 
755 
756 
758 
759 }; // Class Geometry
760 
762 
765 
766 
770 
771 
773 template<class TPointType>
774 inline std::istream& operator >> (std::istream& rIStream,
775  Point3D<TPointType>& rThis);
776 
778 template<class TPointType>
779 inline std::ostream& operator << (std::ostream& rOStream,
780  const Point3D<TPointType>& rThis)
781 {
782  rThis.PrintInfo(rOStream);
783  rOStream << std::endl;
784  rThis.PrintData(rOStream);
785 
786  return rOStream;
787 }
789 
790 template<class TPointType>
791 const GeometryData Point3D<TPointType>::msGeometryData(
794  Point3D<TPointType>::AllIntegrationPoints(),
795  Point3D<TPointType>::AllShapeFunctionsValues(),
796  AllShapeFunctionsLocalGradients());
797 
798 template<class TPointType>
799 const GeometryDimension Point3D<TPointType>::msGeometryDimension(3, 0);
800 
801 } // namespace Kratos.
802 
803 #endif // KRATOS_LINE_2D_H_INCLUDED defined
Definition: geometry_data.h:60
KratosGeometryType
Definition: geometry_data.h:110
IntegrationMethod
Definition: geometry_data.h:76
KratosGeometryFamily
Definition: geometry_data.h:91
Definition: geometry_dimension.h:42
Geometry base class.
Definition: geometry.h:71
SizeType PointsNumber() const
Definition: geometry.h:528
Geometry & operator=(const Geometry &rOther)
Definition: geometry.h:400
std::vector< IntegrationPointType > IntegrationPointsArrayType
Definition: geometry.h:161
DataValueContainer & GetData()
Definition: geometry.h:591
std::size_t SizeType
Definition: geometry.h:144
std::size_t IndexType
Definition: geometry.h:137
std::array< Matrix, static_cast< int >GeometryData::IntegrationMethod::NumberOfIntegrationMethods)> ShapeFunctionsValuesContainerType
Definition: geometry.h:172
const PointsArrayType & Points() const
Definition: geometry.h:1768
std::array< IntegrationPointsArrayType, static_cast< int >GeometryData::IntegrationMethod::NumberOfIntegrationMethods)> IntegrationPointsContainerType
Definition: geometry.h:167
GeometryData::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: geometry.h:177
Short class definition.
Definition: integration_point.h:52
Definition: amatrix_interface.h:41
Definition: point_3d.h:53
Geometry< TPointType > BaseType
Geometry as base class.
Definition: point_3d.h:60
BaseType::ShapeFunctionsValuesContainerType ShapeFunctionsValuesContainerType
Definition: point_3d.h:118
BaseType::CoordinatesArrayType CoordinatesArrayType
Definition: point_3d.h:144
BaseType::IntegrationPointsContainerType IntegrationPointsContainerType
Definition: point_3d.h:113
BaseType::GeometriesArrayType GeometriesArrayType
Definition: point_3d.h:72
BaseType::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: point_3d.h:123
GeometryData::IntegrationMethod IntegrationMethod
Definition: point_3d.h:67
Point3D & operator=(const Point3D &rOther)
Definition: point_3d.h:243
BaseType::Pointer Create(const IndexType NewGeometryId, PointsArrayType const &rThisPoints) const override
Creates a new geometry pointer.
Definition: point_3d.h:278
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
Definition: point_3d.h:220
~Point3D() override
Destructor. Do nothing!!!
Definition: point_3d.h:218
friend class Point3D
Definition: point_3d.h:729
GeometryData::KratosGeometryType GetGeometryType() const override
Definition: point_3d.h:224
Point3D & operator=(Point3D< TOtherPointType > const &rOther)
Definition: point_3d.h:261
SizeType FacesNumber() const override
Returns the number of faces of the current geometry.
Definition: point_3d.h:542
std::string Info() const override
Definition: point_3d.h:594
double Length() const override
Definition: point_3d.h:333
KRATOS_CLASS_POINTER_DEFINITION(Point3D)
Pointer definition of Point3D.
void PrintInfo(std::ostream &rOStream) const override
Definition: point_3d.h:605
BaseType::IntegrationPointType IntegrationPointType
Definition: point_3d.h:100
Point3D(const PointsArrayType &ThisPoints)
Definition: point_3d.h:162
Point3D(const IndexType GeometryId, const PointsArrayType &rThisPoints)
Constructor with Geometry Id.
Definition: point_3d.h:170
double Area() const override
Definition: point_3d.h:349
void PrintData(std::ostream &rOStream) const override
Definition: point_3d.h:618
Point3D(Point3D const &rOther)
Definition: point_3d.h:195
BaseType::Pointer Create(const IndexType NewGeometryId, const BaseType &rGeometry) const override
Creates a new geometry pointer.
Definition: point_3d.h:292
Point3D(Point3D< TOtherPointType > const &rOther)
Definition: point_3d.h:212
BaseType::JacobiansType JacobiansType
Definition: point_3d.h:129
BaseType::IndexType IndexType
Definition: point_3d.h:82
SizeType EdgesNumber() const override
Definition: point_3d.h:537
TPointType PointType
Definition: point_3d.h:76
Point3D(const std::string &rGeometryName, const PointsArrayType &rThisPoints)
Constructor with Geometry Name.
Definition: point_3d.h:179
BaseType::NormalType NormalType
Definition: point_3d.h:139
double DomainSize() const override
Definition: point_3d.h:365
Point3D(typename PointType::Pointer pFirstPoint)
Definition: point_3d.h:156
BaseType::PointsArrayType PointsArrayType
Definition: point_3d.h:94
BaseType::SizeType SizeType
Definition: point_3d.h:89
BaseType::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: point_3d.h:135
BaseType::IntegrationPointsArrayType IntegrationPointsArrayType
Definition: point_3d.h:107
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
void push_back(const TPointerType &x)
Definition: pointer_vector.h:270
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
Short class definition.
Definition: array_1d.h:61
#define KRATOS_SERIALIZE_SAVE_BASE_CLASS(Serializer, BaseType)
Definition: define.h:812
#define KRATOS_SERIALIZE_LOAD_BASE_CLASS(Serializer, BaseType)
Definition: define.h:815
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
const GeometryData Point3D< TPointType >::msGeometryData & msGeometryDimension(), Point3D< TPointType >::AllShapeFunctionsValues(), AllShapeFunctionsLocalGradients()
Definition: brep_curve.h:483
std::istream & operator>>(std::istream &rIStream, LinearMasterSlaveConstraint &rThis)
input stream function
const GeometryData BrepCurve< TContainerPointType, TContainerPointEmbeddedType >::msGeometryData & msGeometryDimension
Definition: brep_curve.h:483
std::ostream & operator<<(std::ostream &rOStream, const LinearMasterSlaveConstraint &rThis)
output stream function
Definition: linear_master_slave_constraint.h:432
def load(f)
Definition: ode_solve.py:307