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.
sphere_3d_1.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_SPHERE_3D_1_H_INCLUDED )
18 #define KRATOS_SPHERE_3D_1_H_INCLUDED
19 
20 // System includes
21 
22 // External includes
23 
24 // Project includes
25 #include "geometries/geometry.h"
27 
28 
29 namespace Kratos
30 {
31 
34 
38 
42 
46 
50 
53 template<class TPointType>
54 
55 class Sphere3D1 : public Geometry<TPointType>
56 {
57 
58 public:
62 
65 
68 
72 
77 
80  typedef TPointType PointType;
81 
86  typedef typename BaseType::IndexType IndexType;
87 
88 
93  typedef typename BaseType::SizeType SizeType;
94 
99 
105 
112 
118 
123 
128 
134 
140 
144 
149 
153 
154 // Sphere3D1( const PointType& FirstPoint)
155 // : BaseType( PointsArrayType(), &msGeometryData )
156 // {
157 // BaseType::Points().push_back( typename PointType::Pointer( new PointType( FirstPoint ) ) );
158 // }
159 
160  Sphere3D1( typename PointType::Pointer pFirstPoint )
161  : BaseType( PointsArrayType(), &msGeometryData )
162  {
163  BaseType::Points().push_back( pFirstPoint );
164  }
165 
166  Sphere3D1( const PointsArrayType& ThisPoints )
167  : BaseType( ThisPoints, &msGeometryData )
168  {
169  if ( BaseType::PointsNumber() != 1 )
170  KRATOS_ERROR << "Invalid points number. Expected 1, given " << BaseType::PointsNumber() << std::endl;
171  }
172 
174  explicit Sphere3D1(
175  const IndexType GeometryId,
176  const PointsArrayType& rThisPoints
177  ) : BaseType(GeometryId, rThisPoints, &msGeometryData)
178  {
179  KRATOS_ERROR_IF( this->PointsNumber() != 1 ) << "Invalid points number. Expected 1, given " << this->PointsNumber() << std::endl;
180  }
181 
183  explicit Sphere3D1(
184  const std::string& rGeometryName,
185  const PointsArrayType& rThisPoints
186  ) : BaseType(rGeometryName, rThisPoints, &msGeometryData)
187  {
188  KRATOS_ERROR_IF(this->PointsNumber() != 1) << "Invalid points number. Expected 1, given " << this->PointsNumber() << std::endl;
189  }
190 
199  Sphere3D1( Sphere3D1 const& rOther )
200  : BaseType( rOther )
201  {
202  }
203 
204 
216  template<class TOtherPointType> Sphere3D1( Sphere3D1<TOtherPointType> const& rOther )
217  : BaseType( rOther )
218  {
219  }
220 
222  ~Sphere3D1() override {}
223 
225  {
227  }
228 
230  {
232  }
233 
237 
248  Sphere3D1& operator=( const Sphere3D1& rOther )
249  {
250  BaseType::operator=( rOther );
251 
252  return *this;
253  }
254 
265  template<class TOtherPointType>
267  {
268  BaseType::operator=( rOther );
269 
270  return *this;
271  }
272 
276 
282  typename BaseType::Pointer Create(
283  PointsArrayType const& rThisPoints
284  ) const override
285  {
286  return typename BaseType::Pointer( new Sphere3D1( rThisPoints ) );
287  }
288 
295  typename BaseType::Pointer Create(
296  const IndexType NewGeometryId,
297  PointsArrayType const& rThisPoints
298  ) const override
299  {
300  return typename BaseType::Pointer( new Sphere3D1( NewGeometryId, rThisPoints ) );
301  }
302 
308  typename BaseType::Pointer Create(
309  const BaseType& rGeometry
310  ) const override
311  {
312  auto p_geometry = typename BaseType::Pointer( new Sphere3D1( rGeometry.Points() ) );
313  p_geometry->SetData(rGeometry.GetData());
314  return p_geometry;
315  }
316 
323  typename BaseType::Pointer Create(
324  const IndexType NewGeometryId,
325  const BaseType& rGeometry
326  ) const override
327  {
328  auto p_geometry = typename BaseType::Pointer( new Sphere3D1( NewGeometryId, rGeometry.Points() ) );
329  p_geometry->SetData(rGeometry.GetData());
330  return p_geometry;
331  }
332 
342  Vector& rResult,
343  const typename BaseType::LumpingMethods LumpingMethod = BaseType::LumpingMethods::ROW_SUM
344  ) const override
345  {
346  if(rResult.size() != 1)
347  rResult.resize( 1, false );
348  rResult[0] = 1.0;
349  return rResult;
350  }
351 
355 
368  double Length() const override
369  {
370  KRATOS_WARNING("Sphere3D1") << "This method (Length) has no meaning for this type of geometry (Sphere)." << std::endl;
371  return 0.0;
372  }
373 
385  double Area() const override
386  {
387  KRATOS_WARNING("Sphere3D1") << "This method (Area) has no meaning for this type of geometry (Sphere)." << std::endl;
388  return 0.0;
389  }
390 
391 
402  double DomainSize() const override
403  {
404  KRATOS_WARNING("Sphere3D1") << "This method (DomainSize) has no meaning for this type of geometry (Sphere)." << std::endl;
405  return 0.0;
406  }
407 
411 
426  JacobiansType& Jacobian( JacobiansType& rResult, IntegrationMethod ThisMethod ) const override
427  {
428  KRATOS_WARNING("Sphere3D1") << "This method (Jacobian) has no meaning for this type of geometry (Sphere)." << std::endl;
429  return rResult;
430  }
431 
449  JacobiansType& Jacobian( JacobiansType& rResult, IntegrationMethod ThisMethod, Matrix & DeltaPosition ) const override
450  {
451  KRATOS_WARNING("Sphere3D1") << "This method (Jacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
452  return rResult;
453  }
454 
472  Matrix& Jacobian( Matrix& rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod ) const override
473  {
474  KRATOS_WARNING("Sphere3D1") << "This method (Jacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
475  return rResult;
476  }
477 
489  Matrix& Jacobian( Matrix& rResult, const CoordinatesArrayType& rPoint ) const override
490  {
491  KRATOS_WARNING("Sphere3D1") << "This method (Jacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
492  return rResult;
493  }
494 
506  Vector& DeterminantOfJacobian( Vector& rResult, IntegrationMethod ThisMethod ) const override
507  {
508  KRATOS_WARNING("Sphere3D1")<<"This method (DeterminantOfJacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
509  return rResult;
510  }
511 
530  double DeterminantOfJacobian( IndexType IntegrationPointIndex, IntegrationMethod ThisMethod ) const override
531  {
532  KRATOS_WARNING("Sphere3D1")<<"This method (DeterminantOfJacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
533  return 0.0;
534  }
535 
548  double DeterminantOfJacobian( const CoordinatesArrayType& rPoint ) const override
549  {
550  KRATOS_WARNING("Sphere3D1")<<"This method (DeterminantOfJacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
551  return 0.0;
552  }
553 
554 
569  JacobiansType& InverseOfJacobian( JacobiansType& rResult, IntegrationMethod ThisMethod ) const override
570  {
571  KRATOS_WARNING("Sphere3D1")<<"This method (InverseOfJacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
572  return rResult;
573  }
574 
592  Matrix& InverseOfJacobian( Matrix& rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod ) const override
593  {
594  KRATOS_WARNING("Sphere3D1") << "This method (InverseOfJacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
595  return rResult;
596  }
597 
609  Matrix& InverseOfJacobian( Matrix& rResult, const CoordinatesArrayType& rPoint ) const override
610  {
611  KRATOS_WARNING("Sphere3D1") << "This method (InverseOfJacobian) has no meaning for this type of geometry (Sphere)."<<std::endl;
612  return rResult;
613  }
614 
615 
619  SizeType EdgesNumber() const override
620  {
621  return 1;
622  }
623 
627 
628  double ShapeFunctionValue( IndexType ShapeFunctionIndex,
629  const CoordinatesArrayType& rPoint ) const override
630  {
631  KRATOS_WARNING("Sphere3D1") << "This method (ShapeFunctionValue) has no meaning for this type of geometry (Sphere)."<<std::endl;
632  return 0;
633  }
634 
636  const CoordinatesArrayType& rPoint ) const override
637  {
638  KRATOS_WARNING("Sphere3D1") << "This method (ShapeFunctionsLocalGradients) has no meaning for this type of geometry (Sphere)."<<std::endl;
639  return( rResult );
640  }
641 
643  {
644  KRATOS_WARNING("Sphere3D1") << "This method (ShapeFunctionsLocalGradients) has no meaning for this type of geometry (Sphere)."<<std::endl;
645  }
646 
653  std::string Info() const override
654  {
655  return "a sphere with 1 nodes in its center, in 3D space";
656  }
657 
664  void PrintInfo( std::ostream& rOStream ) const override
665  {
666  rOStream << "a sphere with 1 nodes in its center, in 3D space";
667  }
668 
677  void PrintData( std::ostream& rOStream ) const override
678  {
679 
680  }
681 
685 
686 
688 
689 protected:
692 
693 
697 
698 
702 
703 
707 
708 
712 
713 
717 
718 
722 
723 
725 
726 private:
729 
730  static const GeometryData msGeometryData;
731 
732  static const GeometryDimension msGeometryDimension;
733 
737 
741 
742  friend class Serializer;
743 
744  void save( Serializer& rSerializer ) const override
745  {
747  }
748 
749  void load( Serializer& rSerializer ) override
750  {
752  }
753 
754  Sphere3D1(): BaseType( PointsArrayType(), &msGeometryData ) {}
755 
756 
760 
761 
765 
766  static Matrix CalculateShapeFunctionsIntegrationPointsValues( typename BaseType::IntegrationMethod ThisMethod )
767  {
768  const IntegrationPointsContainerType& all_integration_points = AllIntegrationPoints();
769  const IntegrationPointsArrayType& IntegrationPoints = all_integration_points[static_cast<int>(ThisMethod)];
770  const int integration_points_number = IntegrationPoints.size();
771  Matrix N( integration_points_number, 1 );
772 
773  //KRATOS_WARNING("Sphere3D1") << "This method (CalculateShapeFunctionsIntegrationPointsValues) has no meaning for this type of geometry (Sphere)."<<std::endl;
774 
775  return N;
776  }
777 
778  static ShapeFunctionsGradientsType CalculateShapeFunctionsIntegrationPointsLocalGradients( typename BaseType::IntegrationMethod ThisMethod )
779  {
780  const IntegrationPointsContainerType& all_integration_points = AllIntegrationPoints();
781  const IntegrationPointsArrayType& IntegrationPoints = all_integration_points[static_cast<int>(ThisMethod)];
783  std::fill( DN_De.begin(), DN_De.end(), Matrix( 2, 1 ) );
784 
785  //KRATOS_WARNING("Sphere3D1") << "This method (CalculateShapeFunctionsIntegrationPointsLocalGradients) has no meaning for this type of geometry (Sphere)."<<std::endl;
786 
787  return DN_De;
788 
789  }
790 
791  static const IntegrationPointsContainerType AllIntegrationPoints()
792  {
793  IntegrationPointsContainerType integration_points = {{
794  Quadrature<LineGaussLegendreIntegrationPoints1, 1, IntegrationPoint<3> >::GenerateIntegrationPoints(),
795  Quadrature<LineGaussLegendreIntegrationPoints2, 1, IntegrationPoint<3> >::GenerateIntegrationPoints(),
796  Quadrature<LineGaussLegendreIntegrationPoints3, 1, IntegrationPoint<3> >::GenerateIntegrationPoints(),
797  Quadrature<LineGaussLegendreIntegrationPoints4, 1, IntegrationPoint<3> >::GenerateIntegrationPoints(),
798  Quadrature<LineGaussLegendreIntegrationPoints5, 1, IntegrationPoint<3> >::GenerateIntegrationPoints()
799  }
800  };
801  return integration_points;
802  }
803 
804  static const ShapeFunctionsValuesContainerType AllShapeFunctionsValues()
805  {
806  ShapeFunctionsValuesContainerType shape_functions_values = {{
807  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsValues( GeometryData::IntegrationMethod::GI_GAUSS_1 ),
808  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsValues( GeometryData::IntegrationMethod::GI_GAUSS_2 ),
809  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsValues( GeometryData::IntegrationMethod::GI_GAUSS_3 ),
810  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsValues( GeometryData::IntegrationMethod::GI_GAUSS_4 ),
811  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsValues( GeometryData::IntegrationMethod::GI_GAUSS_5 )
812  }
813  };
814  return shape_functions_values;
815  }
816 
817  static const ShapeFunctionsLocalGradientsContainerType AllShapeFunctionsLocalGradients()
818  {
819  ShapeFunctionsLocalGradientsContainerType shape_functions_local_gradients = {{
820  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients( GeometryData::IntegrationMethod::GI_GAUSS_1 ),
821  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients( GeometryData::IntegrationMethod::GI_GAUSS_2 ),
822  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients( GeometryData::IntegrationMethod::GI_GAUSS_3 ),
823  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients( GeometryData::IntegrationMethod::GI_GAUSS_4 ),
824  Sphere3D1<TPointType>::CalculateShapeFunctionsIntegrationPointsLocalGradients( GeometryData::IntegrationMethod::GI_GAUSS_5 ),
825 
826  }
827  };
828  return shape_functions_local_gradients;
829  }
830 
834 
835 
839 
840 
844 
845  template<class TOtherPointType> friend class Sphere3D1;
846 
850 
852 
853 }; // Class Geometry
854 
856 
859 
860 
864 
865 
867 template<class TPointType>
868 inline std::istream& operator >> ( std::istream& rIStream,
869  Sphere3D1<TPointType>& rThis );
870 
872 template<class TPointType>
873 inline std::ostream& operator << ( std::ostream& rOStream,
874  const Sphere3D1<TPointType>& rThis )
875 {
876  rThis.PrintInfo( rOStream );
877  rOStream << std::endl;
878  rThis.PrintData( rOStream );
879 
880  return rOStream;
881 }
882 
883 
884 
885 template<class TPointType>
886 const GeometryData Sphere3D1<TPointType>::msGeometryData(
889  Sphere3D1<TPointType>::AllIntegrationPoints(),
890  Sphere3D1<TPointType>::AllShapeFunctionsValues(),
891  AllShapeFunctionsLocalGradients()
892 );
893 
894 template<class TPointType>
895 const GeometryDimension Sphere3D1<TPointType>::msGeometryDimension(3, 1);
896 
897 } // namespace Kratos.
898 
899 #endif // KRATOS_SPHERE_3D_1_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
LumpingMethods
This defines the different methods to compute the lumping methods.
Definition: geometry.h:109
const IntegrationPointsArrayType & IntegrationPoints() const
Definition: geometry.h:2284
GeometryData::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: geometry.h:177
GeometryData::IntegrationMethod IntegrationMethod
Definition: geometry.h:122
Short class definition.
Definition: integration_point.h:52
Definition: amatrix_interface.h:41
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
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
Definition: sphere_3d_1.h:56
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
Definition: sphere_3d_1.h:224
Vector & LumpingFactors(Vector &rResult, const typename BaseType::LumpingMethods LumpingMethod=BaseType::LumpingMethods::ROW_SUM) const override
Lumping factors for the calculation of the lumped mass matrix.
Definition: sphere_3d_1.h:341
BaseType::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: sphere_3d_1.h:139
double Area() const override
Definition: sphere_3d_1.h:385
double Length() const override
Definition: sphere_3d_1.h:368
Matrix & InverseOfJacobian(Matrix &rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const override
Definition: sphere_3d_1.h:592
SizeType EdgesNumber() const override
Definition: sphere_3d_1.h:619
BaseType::SizeType SizeType
Definition: sphere_3d_1.h:93
Matrix & InverseOfJacobian(Matrix &rResult, const CoordinatesArrayType &rPoint) const override
Definition: sphere_3d_1.h:609
BaseType::Pointer Create(const BaseType &rGeometry) const override
Creates a new geometry pointer.
Definition: sphere_3d_1.h:308
Vector & DeterminantOfJacobian(Vector &rResult, IntegrationMethod ThisMethod) const override
Definition: sphere_3d_1.h:506
JacobiansType & InverseOfJacobian(JacobiansType &rResult, IntegrationMethod ThisMethod) const override
Definition: sphere_3d_1.h:569
BaseType::JacobiansType JacobiansType
Definition: sphere_3d_1.h:133
BaseType::PointsArrayType PointsArrayType
Definition: sphere_3d_1.h:98
Sphere3D1(Sphere3D1 const &rOther)
Definition: sphere_3d_1.h:199
~Sphere3D1() override
Destructor. Do nothing!!!
Definition: sphere_3d_1.h:222
void PrintData(std::ostream &rOStream) const override
Definition: sphere_3d_1.h:677
BaseType::Pointer Create(const IndexType NewGeometryId, PointsArrayType const &rThisPoints) const override
It creates a new geometry pointer.
Definition: sphere_3d_1.h:295
Matrix & ShapeFunctionsLocalGradients(Matrix &rResult, const CoordinatesArrayType &rPoint) const override
Definition: sphere_3d_1.h:635
JacobiansType & Jacobian(JacobiansType &rResult, IntegrationMethod ThisMethod, Matrix &DeltaPosition) const override
Definition: sphere_3d_1.h:449
GeometryData::IntegrationMethod IntegrationMethod
Definition: sphere_3d_1.h:71
BaseType::Pointer Create(PointsArrayType const &rThisPoints) const override
It creates a new geometry pointer.
Definition: sphere_3d_1.h:282
Geometry< TPointType > BaseType
Geometry as base class.
Definition: sphere_3d_1.h:64
TPointType PointType
Definition: sphere_3d_1.h:80
double DeterminantOfJacobian(const CoordinatesArrayType &rPoint) const override
Definition: sphere_3d_1.h:548
BaseType::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: sphere_3d_1.h:127
Sphere3D1(const IndexType GeometryId, const PointsArrayType &rThisPoints)
Constructor with Geometry Id.
Definition: sphere_3d_1.h:174
BaseType::NormalType NormalType
Definition: sphere_3d_1.h:143
Matrix & Jacobian(Matrix &rResult, IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const override
Definition: sphere_3d_1.h:472
BaseType::IntegrationPointsArrayType IntegrationPointsArrayType
Definition: sphere_3d_1.h:111
BaseType::CoordinatesArrayType CoordinatesArrayType
Definition: sphere_3d_1.h:148
BaseType::IndexType IndexType
Definition: sphere_3d_1.h:86
BaseType::GeometriesArrayType GeometriesArrayType
Definition: sphere_3d_1.h:76
Sphere3D1(Sphere3D1< TOtherPointType > const &rOther)
Definition: sphere_3d_1.h:216
std::string Info() const override
Definition: sphere_3d_1.h:653
Sphere3D1(const std::string &rGeometryName, const PointsArrayType &rThisPoints)
Constructor with Geometry Name.
Definition: sphere_3d_1.h:183
KRATOS_CLASS_POINTER_DEFINITION(Sphere3D1)
Pointer definition of Sphere3D1.
JacobiansType & Jacobian(JacobiansType &rResult, IntegrationMethod ThisMethod) const override
Definition: sphere_3d_1.h:426
GeometryData::KratosGeometryType GetGeometryType() const override
Definition: sphere_3d_1.h:229
Matrix & Jacobian(Matrix &rResult, const CoordinatesArrayType &rPoint) const override
Definition: sphere_3d_1.h:489
void ShapeFunctionsIntegrationPointsGradients(ShapeFunctionsGradientsType &rResult, IntegrationMethod ThisMethod) const override
Definition: sphere_3d_1.h:642
Sphere3D1(typename PointType::Pointer pFirstPoint)
Definition: sphere_3d_1.h:160
void PrintInfo(std::ostream &rOStream) const override
Definition: sphere_3d_1.h:664
Sphere3D1 & operator=(const Sphere3D1 &rOther)
Definition: sphere_3d_1.h:248
BaseType::IntegrationPointsContainerType IntegrationPointsContainerType
Definition: sphere_3d_1.h:117
Sphere3D1(const PointsArrayType &ThisPoints)
Definition: sphere_3d_1.h:166
double DomainSize() const override
Definition: sphere_3d_1.h:402
double ShapeFunctionValue(IndexType ShapeFunctionIndex, const CoordinatesArrayType &rPoint) const override
Definition: sphere_3d_1.h:628
friend class Sphere3D1
Definition: sphere_3d_1.h:845
Sphere3D1 & operator=(Sphere3D1< TOtherPointType > const &rOther)
Definition: sphere_3d_1.h:266
BaseType::ShapeFunctionsValuesContainerType ShapeFunctionsValuesContainerType
Definition: sphere_3d_1.h:122
BaseType::Pointer Create(const IndexType NewGeometryId, const BaseType &rGeometry) const override
Creates a new geometry pointer.
Definition: sphere_3d_1.h:323
BaseType::IntegrationPointType IntegrationPointType
Definition: sphere_3d_1.h:104
double DeterminantOfJacobian(IndexType IntegrationPointIndex, IntegrationMethod ThisMethod) const override
Definition: sphere_3d_1.h:530
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
const GeometryData Sphere3D1< TPointType >::msGeometryData & msGeometryDimension(), Sphere3D1< TPointType >::AllShapeFunctionsValues(), AllShapeFunctionsLocalGradients()
Definition: brep_curve.h:483
#define KRATOS_ERROR
Definition: exception.h:161
#define KRATOS_ERROR_IF(conditional)
Definition: exception.h:162
#define KRATOS_WARNING(label)
Definition: logger.h:265
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
Internals::Matrix< double, AMatrix::dynamic, AMatrix::dynamic > Matrix
Definition: amatrix_interface.h:470
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
N
Definition: sensitivityMatrix.py:29