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.
geometry_data.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
11 //
12 //
13 
14 #if !defined(KRATOS_GEOMETRY_DATA_H_INCLUDED )
15 #define KRATOS_GEOMETRY_DATA_H_INCLUDED
16 
17 // System includes
18 
19 // External includes
20 
21 // Project includes
26 
27 namespace Kratos
28 {
29 
32 
36 
40 
44 
48 
60 {
61 public:
64 
76  enum class IntegrationMethod {
77  GI_GAUSS_1,
78  GI_GAUSS_2,
79  GI_GAUSS_3,
80  GI_GAUSS_4,
81  GI_GAUSS_5,
87  NumberOfIntegrationMethods // Note that this entry needs to be always the last to be used as integration methods counter
88  };
89 
91  {
101  Kratos_Nurbs,
102  Kratos_Brep,
106  NumberOfGeometryFamilies // Note that this entry needs to be always the last to be used as geometry families counter
107  };
108 
110  {
154  NumberOfGeometryTypes // Note that this entry needs to be always the last to be used as geometry types counter
155  };
159 
162 
167  typedef std::size_t IndexType;
168 
173  typedef std::size_t SizeType;
174 
180 
187 
193 
198 
203 
209 
211 
217 
221 
285  GeometryData(GeometryDimension const *pThisGeometryDimension,
286  IntegrationMethod ThisDefaultMethod,
287  const IntegrationPointsContainerType& ThisIntegrationPoints,
288  const ShapeFunctionsValuesContainerType& ThisShapeFunctionsValues,
289  const ShapeFunctionsLocalGradientsContainerType& ThisShapeFunctionsLocalGradients)
290  : mpGeometryDimension(pThisGeometryDimension)
291  , mGeometryShapeFunctionContainer(
292  ThisDefaultMethod,
293  ThisIntegrationPoints,
294  ThisShapeFunctionsValues,
295  ThisShapeFunctionsLocalGradients)
296  {
297  }
298 
299  /*
300  * Constructor which has a precomputed shape function container.
301  * @param pThisGeometryDimension pointer to the dimensional data
302  * @param ThisGeometryShapeFunctionContainer including the evaluated
303  * values for the shape functions, it's derivatives and the
304  * integration points.
305  */
306  GeometryData(GeometryDimension const *pThisGeometryDimension,
307  const GeometryShapeFunctionContainer<IntegrationMethod>& ThisGeometryShapeFunctionContainer)
308  : mpGeometryDimension(pThisGeometryDimension)
309  , mGeometryShapeFunctionContainer(
310  ThisGeometryShapeFunctionContainer)
311  {
312  }
313 
315  GeometryData( const GeometryData& rOther )
316  : mpGeometryDimension( rOther.mpGeometryDimension)
317  , mGeometryShapeFunctionContainer( rOther.mGeometryShapeFunctionContainer)
318  {
319  }
320 
322  virtual ~GeometryData() {}
323 
327 
330  {
331  mpGeometryDimension = rOther.mpGeometryDimension;
332  mGeometryShapeFunctionContainer = rOther.mGeometryShapeFunctionContainer;
333 
334  return *this;
335  }
336 
340 
341  void SetGeometryDimension(GeometryDimension const* pGeometryDimension)
342  {
343  mpGeometryDimension = pGeometryDimension;
344  }
345 
349 
352  const GeometryShapeFunctionContainer<IntegrationMethod>& rGeometryShapeFunctionContainer)
353  {
354  mGeometryShapeFunctionContainer = rGeometryShapeFunctionContainer;
355  }
356 
359  {
360  return mGeometryShapeFunctionContainer;
361  }
362 
366 
375  {
376  return mpGeometryDimension->WorkingSpaceDimension();
377  }
378 
388  {
389  return mpGeometryDimension->LocalSpaceDimension();
390  }
391 
392 
396 
407  bool HasIntegrationMethod( IntegrationMethod ThisMethod ) const
408  {
409  return mGeometryShapeFunctionContainer.HasIntegrationMethod(ThisMethod);
410  }
411 
415 
426  {
427  return mGeometryShapeFunctionContainer.DefaultIntegrationMethod();
428  }
429 
431  {
432  return mGeometryShapeFunctionContainer.IntegrationPointsNumber();
433  }
434 
444  {
445  return mGeometryShapeFunctionContainer.IntegrationPointsNumber(ThisMethod);
446  }
447 
448 
458  {
459  return mGeometryShapeFunctionContainer.IntegrationPoints();
460  }
461 
471  {
472  return mGeometryShapeFunctionContainer.IntegrationPoints(ThisMethod);
473  }
474 
478 
498  {
499  return mGeometryShapeFunctionContainer.ShapeFunctionsValues();
500  }
501 
522  const Matrix& ShapeFunctionsValues( IntegrationMethod ThisMethod ) const
523  {
524  return mGeometryShapeFunctionContainer.ShapeFunctionsValues(
525  ThisMethod);
526  }
527 
550  double ShapeFunctionValue( IndexType IntegrationPointIndex, IndexType ShapeFunctionIndex ) const
551  {
552  return mGeometryShapeFunctionContainer.ShapeFunctionValue(
553  IntegrationPointIndex,
554  ShapeFunctionIndex);
555  }
556 
580  IndexType IntegrationPointIndex,
581  IndexType ShapeFunctionIndex,
582  IntegrationMethod ThisMethod ) const
583  {
584  return mGeometryShapeFunctionContainer.ShapeFunctionValue(
585  IntegrationPointIndex, ShapeFunctionIndex, ThisMethod );
586  }
587 
608  {
609  return mGeometryShapeFunctionContainer.ShapeFunctionsLocalGradients();
610  }
611 
634  {
635  return mGeometryShapeFunctionContainer.ShapeFunctionsLocalGradients(ThisMethod);
636  }
637 
660  const Matrix& ShapeFunctionLocalGradient( IndexType IntegrationPointIndex ) const
661  {
662  return mGeometryShapeFunctionContainer.ShapeFunctionLocalGradient(IntegrationPointIndex);
663  }
664 
688  const Matrix& ShapeFunctionLocalGradient( IndexType IntegrationPointIndex, IntegrationMethod ThisMethod ) const
689  {
690  return mGeometryShapeFunctionContainer.ShapeFunctionLocalGradient(IntegrationPointIndex, ThisMethod);
691  }
692 
693  const Matrix& ShapeFunctionLocalGradient( IndexType IntegrationPointIndex, IndexType ShapeFunctionIndex, IntegrationMethod ThisMethod ) const
694  {
695  return mGeometryShapeFunctionContainer.ShapeFunctionLocalGradient(IntegrationPointIndex, ThisMethod);
696  }
697 
698  /*
699  * @brief access to the shape function derivatives.
700  * @param DerivativeOrderIndex defines the wanted order of the derivative
701  * @param IntegrationPointIndex the corresponding contorl point of this geometry
702  * @return the shape function or derivative value related to the input parameters
703  * the matrix is structured: (derivative dN_de / dN_du , the corresponding node)
704  */
706  IndexType DerivativeOrderIndex,
707  IndexType IntegrationPointIndex,
708  IntegrationMethod ThisMethod) const
709  {
710  return mGeometryShapeFunctionContainer.ShapeFunctionDerivatives(
711  DerivativeOrderIndex, IntegrationPointIndex, ThisMethod);
712  }
713 
717 
719  virtual std::string Info() const
720  {
721  return "geometry data";
722  }
723 
725  virtual void PrintInfo( std::ostream& rOStream ) const
726  {
727  rOStream << "geometry data";
728  }
729 
731  virtual void PrintData( std::ostream& rOStream ) const
732  {
733  rOStream << " Working space dimension : " << mpGeometryDimension->WorkingSpaceDimension() << std::endl;
734  rOStream << " Local space dimension : " << mpGeometryDimension->LocalSpaceDimension();
735  }
736 
737 
739 
740 private:
743 
744  GeometryDimension const* mpGeometryDimension;
745 
746  GeometryShapeFunctionContainer<IntegrationMethod> mGeometryShapeFunctionContainer;
747 
751 
752  friend class Serializer;
753 
754  virtual void save( Serializer& rSerializer ) const
755  {
756  rSerializer.save("GeometryDimension", mpGeometryDimension);
757  rSerializer.save("GeometryShapeFunctionContainer", mGeometryShapeFunctionContainer);
758  }
759 
760  virtual void load( Serializer& rSerializer )
761  {
762  rSerializer.load("GeometryDimension", const_cast<GeometryDimension*>(mpGeometryDimension));
763  rSerializer.load("GeometryShapeFunctionContainer", mGeometryShapeFunctionContainer);
764  }
765 
766  // Private default constructor for serialization
767  GeometryData()
768  {
769  }
770 
772 
773 }; // Class GeometryData
774 
776 
779 
780 
784 
785 
787 inline std::istream& operator >> ( std::istream& rIStream,
788  GeometryData& rThis );
789 
791 inline std::ostream& operator << ( std::ostream& rOStream,
792  const GeometryData& rThis )
793 {
794  rThis.PrintInfo( rOStream );
795  rOStream << std::endl;
796  rThis.PrintData( rOStream );
797 
798  return rOStream;
799 }
800 
802 
803 } // namespace Kratos.
804 
805 #endif // KRATOS_GEOMETRY_DATA_H_INCLUDED defined
Definition: geometry_data.h:60
GeometryShapeFunctionContainer< IntegrationMethod >::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: geometry_data.h:202
const ShapeFunctionsGradientsType & ShapeFunctionsLocalGradients() const
Definition: geometry_data.h:607
double ShapeFunctionValue(IndexType IntegrationPointIndex, IndexType ShapeFunctionIndex) const
Definition: geometry_data.h:550
IntegrationPoint< 3 > IntegrationPointType
Definition: geometry_data.h:179
const Matrix & ShapeFunctionLocalGradient(IndexType IntegrationPointIndex) const
Definition: geometry_data.h:660
const Matrix & ShapeFunctionLocalGradient(IndexType IntegrationPointIndex, IndexType ShapeFunctionIndex, IntegrationMethod ThisMethod) const
Definition: geometry_data.h:693
KratosGeometryType
Definition: geometry_data.h:110
SizeType IntegrationPointsNumber(IntegrationMethod ThisMethod) const
Definition: geometry_data.h:443
GeometryShapeFunctionContainer< IntegrationMethod >::IntegrationPointsContainerType IntegrationPointsContainerType
Definition: geometry_data.h:192
const ShapeFunctionsGradientsType & ShapeFunctionsLocalGradients(IntegrationMethod ThisMethod) const
Definition: geometry_data.h:633
SizeType LocalSpaceDimension() const
Definition: geometry_data.h:387
void SetGeometryShapeFunctionContainer(const GeometryShapeFunctionContainer< IntegrationMethod > &rGeometryShapeFunctionContainer)
SetGeometryShapeFunctionContainer updates the GeometryShapeFunctionContainer.
Definition: geometry_data.h:351
SizeType IntegrationPointsNumber() const
Definition: geometry_data.h:430
GeometryShapeFunctionContainer< IntegrationMethod >::ShapeFunctionsValuesContainerType ShapeFunctionsValuesContainerType
Definition: geometry_data.h:197
virtual void PrintInfo(std::ostream &rOStream) const
Print information about this object.
Definition: geometry_data.h:725
IntegrationMethod
Definition: geometry_data.h:76
bool HasIntegrationMethod(IntegrationMethod ThisMethod) const
Definition: geometry_data.h:407
const Matrix & ShapeFunctionLocalGradient(IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const
Definition: geometry_data.h:688
SizeType WorkingSpaceDimension() const
Definition: geometry_data.h:374
const GeometryShapeFunctionContainer< IntegrationMethod > & GetGeometryShapeFunctionContainer() const
Returns the GeometryShapeFunctionContainer.
Definition: geometry_data.h:358
DenseVector< DenseVector< Matrix > > ShapeFunctionsThirdDerivativesType
Definition: geometry_data.h:216
const Matrix & ShapeFunctionsValues() const
Definition: geometry_data.h:497
void SetGeometryDimension(GeometryDimension const *pGeometryDimension)
Definition: geometry_data.h:341
DenseVector< Matrix > ShapeFunctionsGradientsType
Definition: geometry_data.h:208
std::size_t SizeType
Definition: geometry_data.h:173
virtual void PrintData(std::ostream &rOStream) const
Print object's data.
Definition: geometry_data.h:731
virtual std::string Info() const
Turn back information as a string.
Definition: geometry_data.h:719
GeometryShapeFunctionContainer< IntegrationMethod >::IntegrationPointsArrayType IntegrationPointsArrayType
Definition: geometry_data.h:186
const Matrix & ShapeFunctionsValues(IntegrationMethod ThisMethod) const
Definition: geometry_data.h:522
GeometryData(GeometryDimension const *pThisGeometryDimension, IntegrationMethod ThisDefaultMethod, const IntegrationPointsContainerType &ThisIntegrationPoints, const ShapeFunctionsValuesContainerType &ThisShapeFunctionsValues, const ShapeFunctionsLocalGradientsContainerType &ThisShapeFunctionsLocalGradients)
Definition: geometry_data.h:285
const IntegrationPointsArrayType & IntegrationPoints() const
Definition: geometry_data.h:457
KRATOS_CLASS_POINTER_DEFINITION(GeometryData)
Pointer definition of GeometryData.
KratosGeometryFamily
Definition: geometry_data.h:91
const IntegrationPointsArrayType & IntegrationPoints(IntegrationMethod ThisMethod) const
Definition: geometry_data.h:470
std::size_t IndexType
Definition: geometry_data.h:167
double ShapeFunctionValue(IndexType IntegrationPointIndex, IndexType ShapeFunctionIndex, IntegrationMethod ThisMethod) const
Definition: geometry_data.h:579
GeometryData & operator=(const GeometryData &rOther)
Assignment operator.
Definition: geometry_data.h:329
const Matrix & ShapeFunctionDerivatives(IndexType DerivativeOrderIndex, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const
Definition: geometry_data.h:705
GeometryData(const GeometryData &rOther)
Copy constructor.
Definition: geometry_data.h:315
virtual ~GeometryData()
Destructor.
Definition: geometry_data.h:322
GeometryData(GeometryDimension const *pThisGeometryDimension, const GeometryShapeFunctionContainer< IntegrationMethod > &ThisGeometryShapeFunctionContainer)
Definition: geometry_data.h:306
DenseVector< Matrix > ShapeFunctionsSecondDerivativesType
Definition: geometry_data.h:210
IntegrationMethod DefaultIntegrationMethod() const
Definition: geometry_data.h:425
Definition: geometry_dimension.h:42
SizeType LocalSpaceDimension() const
Definition: geometry_dimension.h:128
SizeType WorkingSpaceDimension() const
Definition: geometry_dimension.h:115
Definition: geometry_shape_function_container.h:60
const Matrix & ShapeFunctionsValues() const
Definition: geometry_shape_function_container.h:245
const IntegrationPointsArrayType & IntegrationPoints() const
Definition: geometry_shape_function_container.h:231
const ShapeFunctionsGradientsType & ShapeFunctionsLocalGradients() const
Definition: geometry_shape_function_container.h:275
SizeType IntegrationPointsNumber() const
Definition: geometry_shape_function_container.h:221
bool HasIntegrationMethod(IntegrationMethod ThisMethod) const
Definition: geometry_shape_function_container.h:212
double ShapeFunctionValue(IndexType IntegrationPointIndex, IndexType ShapeFunctionIndex, IntegrationMethod ThisMethod) const
Definition: geometry_shape_function_container.h:256
IntegrationMethod DefaultIntegrationMethod() const
Definition: geometry_shape_function_container.h:207
const Matrix & ShapeFunctionLocalGradient(IndexType IntegrationPointIndex) const
Definition: geometry_shape_function_container.h:286
const Matrix & ShapeFunctionDerivatives(IndexType DerivativeOrderIndex, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const
Definition: geometry_shape_function_container.h:330
Short class definition.
Definition: integration_point.h:52
Definition: amatrix_interface.h:41
The serialization consists in storing the state of an object into a storage format like data file or ...
Definition: serializer.h:123
void load(std::string const &rTag, TDataType &rObject)
Definition: serializer.h:207
void save(std::string const &rTag, std::array< TDataType, TDataSize > const &rObject)
Definition: serializer.h:545
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
def load(f)
Definition: ode_solve.py:307