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.
quadrature_point_surface_in_volume_geometry.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: Manuel Messmer
11 //
12 
13 #if !defined(KRATOS_QUADRATURE_POINT_SURFACE_IN_VOLUME_GEOMETRY_H_INCLUDED )
14 #define KRATOS_QUADRATURE_POINT_SURFACE_IN_VOLUME_GEOMETRY_H_INCLUDED
15 
16 // System includes
17 
18 // External includes
19 
20 // Project includes
21 #include "includes/variables.h"
23 
24 namespace Kratos
25 {
37 template<class TPointType>
39  : public QuadraturePointGeometry<TPointType, 3, 3, 2>
40 {
41 public:
42 
45 
48 
50  typedef typename GeometryType::SizeType SizeType;
51 
54 
57 
59  using BaseType::Jacobian;
60  using BaseType::Calculate;
62 
66 
69  const PointsArrayType& ThisPoints,
70  GeometryShapeFunctionContainerType& ThisGeometryShapeFunctionContainer,
71  const TangentMatrixType& LocalTangents )
72  : BaseType(ThisPoints, ThisGeometryShapeFunctionContainer)
73  {
74  KRATOS_ERROR_IF( mLocalTangents.size1() != LocalTangents.size1() || mLocalTangents.size2() != LocalTangents.size2() )
75  << "QuadraturePointSurfaceInVolumeGeometry :: Dimensions of LocalTangents do not match (3,2). "
76  << "Given Dimensions are: (" << LocalTangents.size1() << "," << LocalTangents.size2() <<"). " << std::endl;
77  mLocalTangents = LocalTangents;
78  }
79 
82  const PointsArrayType& ThisPoints,
83  GeometryShapeFunctionContainerType& ThisGeometryShapeFunctionContainer,
84  const TangentMatrixType& LocalTangents,
85  GeometryType* pGeometryParent)
86  : BaseType(ThisPoints, ThisGeometryShapeFunctionContainer, pGeometryParent)
87  {
88  KRATOS_ERROR_IF( mLocalTangents.size1() != LocalTangents.size1() || mLocalTangents.size2() != LocalTangents.size2() )
89  << "QuadraturePointSurfaceInVolumeGeometry :: Dimensions of LocalTangents do not match (3,2). "
90  << "Given Dimensions are: (" << LocalTangents.size1() << "," << LocalTangents.size2() <<"). " << std::endl;
91  mLocalTangents = LocalTangents;
92  }
93 
96 
99  : BaseType( rOther )
100  , mLocalTangents(rOther.mLocalTangents)
101  {
102  }
103 
107 
110  {
111  BaseType::operator=( rOther );
112 
113  mLocalTangents = rOther.mLocalTangents;
114 
115  return *this;
116  }
117 
121 
123  void Assign(
124  const Variable<Matrix >& rVariable,
125  const Matrix& rInput) override
126  {
127  KRATOS_ERROR_IF( rInput.size1() != mLocalTangents.size1() || rInput.size2() != mLocalTangents.size2() )
128  << "QuadraturePointSurfaceInVolume :: Input Matrix does not have the same size as LocalTangentMatrix. Size must be (3,2)." << std::endl;
129 
130  if (rVariable == LOCAL_TANGENT_MATRIX) {
131  mLocalTangents = rInput;
132  }
133  }
134 
136  void Calculate(
137  const Variable<Matrix >& rVariable,
138  Matrix& rOutput) const override
139  {
140  if (rVariable == LOCAL_TANGENT_MATRIX) {
141  if( rOutput.size1() != mLocalTangents.size1() || rOutput.size2() != mLocalTangents.size2() ){
142  rOutput.resize(mLocalTangents.size1(),mLocalTangents.size2());
143  }
144  rOutput = mLocalTangents;
145  }
146  }
147 
151 
159  IndexType IntegrationPointIndex,
160  GeometryData::IntegrationMethod ThisMethod) const override
161  {
162  KRATOS_DEBUG_ERROR_IF(IntegrationPointIndex != 0)
163  << "Trying to access Normal of QuadraturePointCurveOnSurface "
164  << "with an integration point index != 0." << std::endl;
165 
166  Matrix J;
167  this->Jacobian(J, IntegrationPointIndex, ThisMethod);
168 
169  // Map tangents into global space
170  TangentMatrixType global_tangents = prod(J, mLocalTangents);
171 
172  Vector g1 = column(global_tangents,0);
173  Vector g2 = column(global_tangents,1);
174 
175  // Compute normal
177 
178  return normal;
179  }
180 
184 
191  IndexType IntegrationPointIndex,
192  GeometryData::IntegrationMethod ThisMethod) const override
193  {
194  KRATOS_DEBUG_ERROR_IF(IntegrationPointIndex != 0)
195  << "Trying to access DeterminantOfJacobian of QuadraturePointSurfaceInVolume "
196  << "with an integration point index != 0." << std::endl;
197 
198  Matrix J;
199  this->Jacobian(J, IntegrationPointIndex, ThisMethod);
200  TangentMatrixType global_tangents = prod(J, mLocalTangents);
201 
202  return MathUtils<double>::GeneralizedDet( global_tangents );
203  }
204 
211  Vector& rResult,
212  GeometryData::IntegrationMethod ThisMethod) const override
213  {
214  if (rResult.size() != 1)
215  rResult.resize(1, false);
216 
217  rResult[0] = this->DeterminantOfJacobian(0, ThisMethod);
218 
219  return rResult;
220  }
221 
225 
227  {
229  }
230 
232  {
234  }
235 
239 
241  std::string Info() const override
242  {
243  std::stringstream buffer;
244  buffer << "Quadrature point geometry for a surface in a nurbs volume with Id: "
245  << std::to_string(this->Id()) << ", containing: " << std::to_string(this->size()) << " points." << std::endl;
246 
247  return buffer.str();
248  }
249 
251  void PrintInfo( std::ostream& rOStream ) const override
252  {
253  rOStream << "Quadrature point geometry for a surface in a nurbs volume with Id: "
254  << this->Id() << ", containing: " << this->size() << " points." << std::endl;
255  }
256 
258  void PrintData( std::ostream& rOStream ) const override
259  {
260  }
262 
263 private:
266 
267  TangentMatrixType mLocalTangents;
268 
272 
275  : BaseType()
276  {
277  }
278 
279  friend class Serializer;
280 
281  void save( Serializer& rSerializer ) const override
282  {
284  rSerializer.save("LocalTangents", mLocalTangents);
285  }
286 
287  void load( Serializer& rSerializer ) override
288  {
290  rSerializer.load("LocalTangents", mLocalTangents);
291  }
292 
294 }; // Class Geometry
295 
296 } // namespace Kratos.
297 
298 #endif // KRATOS_QUADRATURE_POINT_SURFACE_IN_VOLUME_GEOMETRY_H_INCLUDED defined
KratosGeometryType
Definition: geometry_data.h:110
IntegrationMethod
Definition: geometry_data.h:76
KratosGeometryFamily
Definition: geometry_data.h:91
Geometry base class.
Definition: geometry.h:71
SizeType size() const
Definition: geometry.h:518
std::size_t SizeType
Definition: geometry.h:144
std::size_t IndexType
Definition: geometry.h:137
IndexType const & Id() const
Id of this Geometry.
Definition: geometry.h:964
Definition: geometry_shape_function_container.h:60
Definition: amatrix_interface.h:41
void resize(std::size_t NewSize1, std::size_t NewSize2, bool preserve=0)
Definition: amatrix_interface.h:224
static T CrossProduct(const T &a, const T &b)
Performs the vector product of the two input vectors a,b.
Definition: math_utils.h:762
static double GeneralizedDet(const TMatrixType &rA)
Calculates the determinant of any matrix (no check performed on matrix size)
Definition: math_utils.h:634
PointerVector is a container like stl vector but using a vector to store pointers to its data.
Definition: pointer_vector.h:72
A sinlge quadrature point, that can be used for geometries without a predefined integration scheme,...
Definition: quadrature_point_geometry.h:44
Vector & DeterminantOfJacobian(Vector &rResult) const
Definition: geometry.h:3154
QuadraturePointGeometry & operator=(const QuadraturePointGeometry &rOther)
Assignment operator.
Definition: quadrature_point_geometry.h:220
void Calculate(const Variable< array_1d< double, 3 >> &rVariable, array_1d< double, 3 > &rOutput) const override
Calculate with array_1d<double, 3>
Definition: quadrature_point_geometry.h:282
JacobiansType & Jacobian(JacobiansType &rResult) const
using base class functions
Definition: geometry.h:2901
A sinlge quadrature point, that can be used for geometries without a predefined integration scheme,...
Definition: quadrature_point_surface_in_volume_geometry.h:40
void Calculate(const Variable< Matrix > &rVariable, Matrix &rOutput) const override
Calculate with Matrix.
Definition: quadrature_point_surface_in_volume_geometry.h:136
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: quadrature_point_surface_in_volume_geometry.h:258
KRATOS_CLASS_POINTER_DEFINITION(QuadraturePointSurfaceInVolumeGeometry)
Pointer definition of QuadraturePointGeometry.
GeometryData::KratosGeometryType GetGeometryType() const override
Definition: quadrature_point_surface_in_volume_geometry.h:231
GeometryType::PointsArrayType PointsArrayType
Definition: quadrature_point_surface_in_volume_geometry.h:52
QuadraturePointGeometry< TPointType, 3, 3, 2 > BaseType
Definition: quadrature_point_surface_in_volume_geometry.h:46
GeometryType::CoordinatesArrayType CoordinatesArrayType
Definition: quadrature_point_surface_in_volume_geometry.h:53
GeometryShapeFunctionContainer< GeometryData::IntegrationMethod > GeometryShapeFunctionContainerType
Definition: quadrature_point_surface_in_volume_geometry.h:55
BoundedMatrix< double, 3, 2 > TangentMatrixType
Definition: quadrature_point_surface_in_volume_geometry.h:56
QuadraturePointSurfaceInVolumeGeometry(const PointsArrayType &ThisPoints, GeometryShapeFunctionContainerType &ThisGeometryShapeFunctionContainer, const TangentMatrixType &LocalTangents, GeometryType *pGeometryParent)
Constructor with points, geometry shape function container and parent.
Definition: quadrature_point_surface_in_volume_geometry.h:81
QuadraturePointSurfaceInVolumeGeometry(QuadraturePointSurfaceInVolumeGeometry const &rOther)
Copy constructor.
Definition: quadrature_point_surface_in_volume_geometry.h:98
CoordinatesArrayType Normal(IndexType IntegrationPointIndex, GeometryData::IntegrationMethod ThisMethod) const override
Computes the normal of the surface lying on the underlying volume.
Definition: quadrature_point_surface_in_volume_geometry.h:158
QuadraturePointSurfaceInVolumeGeometry & operator=(const QuadraturePointSurfaceInVolumeGeometry &rOther)
Assignment operator.
Definition: quadrature_point_surface_in_volume_geometry.h:109
std::string Info() const override
Turn back information as a string.
Definition: quadrature_point_surface_in_volume_geometry.h:241
Geometry< TPointType > GeometryType
Definition: quadrature_point_surface_in_volume_geometry.h:47
GeometryType::SizeType SizeType
Definition: quadrature_point_surface_in_volume_geometry.h:50
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
Definition: quadrature_point_surface_in_volume_geometry.h:226
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: quadrature_point_surface_in_volume_geometry.h:251
GeometryType::IndexType IndexType
Definition: quadrature_point_surface_in_volume_geometry.h:49
~QuadraturePointSurfaceInVolumeGeometry() override=default
Destructor.
Vector & DeterminantOfJacobian(Vector &rResult, GeometryData::IntegrationMethod ThisMethod) const override
Returns the determinant of the jacobian at this quadrature point.
Definition: quadrature_point_surface_in_volume_geometry.h:210
void Assign(const Variable< Matrix > &rVariable, const Matrix &rInput) override
Assign with Matrix.
Definition: quadrature_point_surface_in_volume_geometry.h:123
double DeterminantOfJacobian(IndexType IntegrationPointIndex, GeometryData::IntegrationMethod ThisMethod) const override
Returns the respective surface area of this quadrature point in global space.
Definition: quadrature_point_surface_in_volume_geometry.h:190
QuadraturePointSurfaceInVolumeGeometry(const PointsArrayType &ThisPoints, GeometryShapeFunctionContainerType &ThisGeometryShapeFunctionContainer, const TangentMatrixType &LocalTangents)
Constructor with points and geometry shape function container.
Definition: quadrature_point_surface_in_volume_geometry.h:68
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
Variable class contains all information needed to store and retrive data from a data container.
Definition: variable.h:63
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_IF(conditional)
Definition: exception.h:162
#define KRATOS_DEBUG_ERROR_IF(conditional)
Definition: exception.h:171
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
AMatrix::MatrixProductExpression< TExpression1Type, TExpression2Type > prod(AMatrix::MatrixExpression< TExpression1Type, TCategory1 > const &First, AMatrix::MatrixExpression< TExpression2Type, TCategory2 > const &Second)
Definition: amatrix_interface.h:568
AMatrix::MatrixColumn< const TExpressionType > column(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression, std::size_t ColumnIndex)
Definition: amatrix_interface.h:637
def load(f)
Definition: ode_solve.py:307
J
Definition: sensitivityMatrix.py:58