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_curve_on_surface_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: Tobias Teschemacher
11 // Veronika Singer
12 // Philipp Bucher
13 //
14 
15 #if !defined(KRATOS_QUADRATURE_POINT_CURVE_ON_SURFACE_GEOMETRY_H_INCLUDED )
16 #define KRATOS_QUADRATURE_POINT_CURVE_ON_SURFACE_GEOMETRY_H_INCLUDED
17 
18 // System includes
19 
20 // External includes
21 
22 // Project includes
23 #include "includes/variables.h"
25 
26 namespace Kratos
27 {
40 template<class TPointType>
42  : public QuadraturePointGeometry<TPointType, 3, 2, 1>
43 {
44 public:
45 
48 
51 
53  typedef typename GeometryType::SizeType SizeType;
54 
57 
60 
62 
64 
68 
70  using BaseType::Jacobian;
75 
79 
82  const PointsArrayType& ThisPoints,
83  GeometryShapeFunctionContainerType& ThisGeometryShapeFunctionContainer,
84  double LocalTangentsU,
85  double LocalTangentsV)
86  : BaseType(ThisPoints, ThisGeometryShapeFunctionContainer)
87  , mLocalTangentsU(LocalTangentsU)
88  , mLocalTangentsV(LocalTangentsV)
89  {
90  }
91 
94  const PointsArrayType& ThisPoints,
95  GeometryShapeFunctionContainerType& ThisGeometryShapeFunctionContainer,
96  double LocalTangentsU,
97  double LocalTangentsV,
98  GeometryType* pGeometryParent)
99  : BaseType(ThisPoints, ThisGeometryShapeFunctionContainer, pGeometryParent)
100  , mLocalTangentsU(LocalTangentsU)
101  , mLocalTangentsV(LocalTangentsV)
102  {
103  }
104 
107 
110  : BaseType( rOther )
111  , mLocalTangentsU(rOther.mLocalTangentsU)
112  , mLocalTangentsV(rOther.mLocalTangentsV)
113  {
114  }
115 
119 
122  {
123  BaseType::operator=( rOther );
124 
125  mLocalTangentsU = rOther.mLocalTangentsU;
126  mLocalTangentsV = rOther.mLocalTangentsV;
127 
128  return *this;
129  }
130 
134 
136  void Assign(
137  const Variable<array_1d<double, 3>>& rVariable,
138  const array_1d<double, 3>& rInput) override
139  {
140  if (rVariable == LOCAL_TANGENT)
141  {
142  mLocalTangentsU = rInput[0];
143  mLocalTangentsV = rInput[1];
144  }
145  }
146 
148  void Calculate(
149  const Variable<array_1d<double, 3>>& rVariable,
150  array_1d<double, 3>& rOutput) const override
151  {
152  if (rVariable == LOCAL_TANGENT)
153  {
154  rOutput[0] = mLocalTangentsU;
155  rOutput[1] = mLocalTangentsV;
156  rOutput[2] = 0.0;
157  }
158  }
159 
161  void Calculate(
162  const Variable<Vector>& rVariable,
163  Vector& rOutput) const override
164  {
165  if (rVariable == DETERMINANTS_OF_JACOBIAN_PARENT)
166  {
168  }
169  }
170 
174 
175  /*
176  * @brief computes the normal of the curve
177  * laying on the underlying surface.
178  * @param rResult Normal to the curve lying on the surface.
179  * @param IntegrationPointIndex index should be always 0.
180  */
182  IndexType IntegrationPointIndex,
183  GeometryData::IntegrationMethod ThisMethod) const override
184  {
185  KRATOS_DEBUG_ERROR_IF(IntegrationPointIndex != 0)
186  << "Trying to access Normal of QuadraturePointCurveOnSurface "
187  << "with an integration point index != 0." << std::endl;
188 
189  Matrix J;
190  this->Jacobian(J, IntegrationPointIndex, ThisMethod);
191 
192  array_1d<double, 3> a_1 = column(J, 0);
193  array_1d<double, 3> a_2 = column(J, 1);
194 
195  CoordinatesArrayType normal = a_2 * mLocalTangentsU - a_1 * mLocalTangentsV;
196 
197  return normal;
198  }
199 
203 
204  /*
205  * @brief returns the respective curve length of this
206  * quadrature point.
207  * @param IntegrationPointIndex index should be always 0.
208  */
210  IndexType IntegrationPointIndex,
211  GeometryData::IntegrationMethod ThisMethod) const override
212  {
213  KRATOS_DEBUG_ERROR_IF(IntegrationPointIndex != 0)
214  << "Trying to access DeterminantOfJacobian of QuadraturePointCurveOnSurface "
215  << "with an integration point index != 0." << std::endl;
216 
217  Matrix J;
218  this->Jacobian(J, IntegrationPointIndex, ThisMethod);
219 
220  array_1d<double, 3> a_1 = column(J, 0);
221  array_1d<double, 3> a_2 = column(J, 1);
222 
223  return norm_2(a_1 * mLocalTangentsU + a_2 * mLocalTangentsV);
224  }
225 
226  /* @brief returns the respective segment length of this
227  * quadrature point. Length of vector always 1.
228  * @param rResult vector of results of this quadrature point.
229  */
231  Vector& rResult,
232  GeometryData::IntegrationMethod ThisMethod) const override
233  {
234  if (rResult.size() != 1)
235  rResult.resize(1, false);
236 
237  rResult[0] = this->DeterminantOfJacobian(0, ThisMethod);
238 
239  return rResult;
240  }
241 
242  /* @brief returns the respective segment length of this
243  * quadrature point, computed on the parent of this geometry.
244  * Required for reduced quadrature point geometries (Not all
245  * nodes are part of this geometry - used for mapping).
246  * @param rResult vector of results of this quadrature point.
247  */
249  Vector& rResult) const
250  {
251  if (rResult.size() != 1)
252  rResult.resize(1, false);
253 
254  Matrix J;
255  this->GetGeometryParent(0).Jacobian(J, this->IntegrationPoints()[0]);
256 
257  rResult[0] = norm_2(column(J, 0) * mLocalTangentsU + column(J, 1) * mLocalTangentsV);
258 
259  return rResult;
260  }
261 
265 
267  {
269  }
270 
272  {
274  }
275 
279 
281  std::string Info() const override
282  {
283  return "Quadrature point for a curve on surface.";
284  }
285 
287  void PrintInfo( std::ostream& rOStream ) const override
288  {
289  rOStream << "Quadrature point for a curve on surface.";
290  }
291 
293  void PrintData( std::ostream& rOStream ) const override
294  {
295  }
297 
298 private:
301 
302  double mLocalTangentsU;
303  double mLocalTangentsV;
304 
308 
311  : BaseType()
312  {
313  }
314 
315  friend class Serializer;
316 
317  void save( Serializer& rSerializer ) const override
318  {
320  rSerializer.save("LocalTangentsU", mLocalTangentsU);
321  rSerializer.save("LocalTangentsV", mLocalTangentsV);
322  }
323 
324  void load( Serializer& rSerializer ) override
325  {
327  rSerializer.load("LocalTangentsU", mLocalTangentsU);
328  rSerializer.load("LocalTangentsV", mLocalTangentsV);
329  }
330 
332 }; // Class Geometry
333 
334 } // namespace Kratos.
335 
336 #endif // KRATOS_QUADRATURE_POINT_CURVE_ON_SURFACE_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
std::vector< IntegrationPointType > IntegrationPointsArrayType
Definition: geometry.h:161
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
std::array< IntegrationPointsArrayType, static_cast< int >GeometryData::IntegrationMethod::NumberOfIntegrationMethods)> IntegrationPointsContainerType
Definition: geometry.h:167
const IntegrationPointsArrayType & IntegrationPoints() const
Definition: geometry.h:2284
JacobiansType & Jacobian(JacobiansType &rResult) const
Definition: geometry.h:2901
GeometryData::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: geometry.h:177
Definition: geometry_shape_function_container.h:60
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
Definition: quadrature_point_curve_on_surface_geometry.h:43
std::string Info() const override
Turn back information as a string.
Definition: quadrature_point_curve_on_surface_geometry.h:281
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_curve_on_surface_geometry.h:148
GeometryData::KratosGeometryFamily GetGeometryFamily() const override
Definition: quadrature_point_curve_on_surface_geometry.h:266
~QuadraturePointCurveOnSurfaceGeometry() override=default
Destructor.
Vector & DeterminantOfJacobianParent(Vector &rResult) const
Definition: quadrature_point_curve_on_surface_geometry.h:248
GeometryType::ShapeFunctionsValuesContainerType ShapeFunctionsValuesContainerType
Definition: quadrature_point_curve_on_surface_geometry.h:66
KRATOS_CLASS_POINTER_DEFINITION(QuadraturePointCurveOnSurfaceGeometry)
Pointer definition of QuadraturePointGeometry.
void PrintInfo(std::ostream &rOStream) const override
Print information about this object.
Definition: quadrature_point_curve_on_surface_geometry.h:287
QuadraturePointCurveOnSurfaceGeometry(const PointsArrayType &ThisPoints, GeometryShapeFunctionContainerType &ThisGeometryShapeFunctionContainer, double LocalTangentsU, double LocalTangentsV)
Constructor with points and geometry shape function container.
Definition: quadrature_point_curve_on_surface_geometry.h:81
CoordinatesArrayType Normal(IndexType IntegrationPointIndex, GeometryData::IntegrationMethod ThisMethod) const override
It returns the vector, which is normal to its corresponding geometry in the given integration point.
Definition: quadrature_point_curve_on_surface_geometry.h:181
GeometryType::IntegrationPointsContainerType IntegrationPointsContainerType
Definition: quadrature_point_curve_on_surface_geometry.h:65
GeometryShapeFunctionContainer< GeometryData::IntegrationMethod > GeometryShapeFunctionContainerType
Definition: quadrature_point_curve_on_surface_geometry.h:63
void Calculate(const Variable< Vector > &rVariable, Vector &rOutput) const override
Calculate with Vector.
Definition: quadrature_point_curve_on_surface_geometry.h:161
GeometryType::IntegrationPointsArrayType IntegrationPointsArrayType
Definition: quadrature_point_curve_on_surface_geometry.h:59
Vector & DeterminantOfJacobian(Vector &rResult, GeometryData::IntegrationMethod ThisMethod) const override
Definition: quadrature_point_curve_on_surface_geometry.h:230
GeometryType::IntegrationPointType IntegrationPointType
Definition: quadrature_point_curve_on_surface_geometry.h:58
void PrintData(std::ostream &rOStream) const override
Print object's data.
Definition: quadrature_point_curve_on_surface_geometry.h:293
GeometryType::IndexType IndexType
Definition: quadrature_point_curve_on_surface_geometry.h:52
GeometryType::SizeType SizeType
Definition: quadrature_point_curve_on_surface_geometry.h:53
GeometryType::CoordinatesArrayType CoordinatesArrayType
Definition: quadrature_point_curve_on_surface_geometry.h:56
QuadraturePointCurveOnSurfaceGeometry & operator=(const QuadraturePointCurveOnSurfaceGeometry &rOther)
Assignment operator.
Definition: quadrature_point_curve_on_surface_geometry.h:121
GeometryData::ShapeFunctionsGradientsType ShapeFunctionsGradientsType
Definition: quadrature_point_curve_on_surface_geometry.h:61
Geometry< TPointType > GeometryType
Definition: quadrature_point_curve_on_surface_geometry.h:50
void Assign(const Variable< array_1d< double, 3 >> &rVariable, const array_1d< double, 3 > &rInput) override
Assign with array_1d<double, 3>
Definition: quadrature_point_curve_on_surface_geometry.h:136
GeometryData::KratosGeometryType GetGeometryType() const override
Definition: quadrature_point_curve_on_surface_geometry.h:271
QuadraturePointCurveOnSurfaceGeometry(const PointsArrayType &ThisPoints, GeometryShapeFunctionContainerType &ThisGeometryShapeFunctionContainer, double LocalTangentsU, double LocalTangentsV, GeometryType *pGeometryParent)
Constructor with points, geometry shape function container, parent.
Definition: quadrature_point_curve_on_surface_geometry.h:93
QuadraturePointGeometry< TPointType, 3, 2, 1 > BaseType
Definition: quadrature_point_curve_on_surface_geometry.h:49
GeometryType::ShapeFunctionsLocalGradientsContainerType ShapeFunctionsLocalGradientsContainerType
Definition: quadrature_point_curve_on_surface_geometry.h:67
double DeterminantOfJacobian(IndexType IntegrationPointIndex, GeometryData::IntegrationMethod ThisMethod) const override
Definition: quadrature_point_curve_on_surface_geometry.h:209
GeometryType::PointsArrayType PointsArrayType
Definition: quadrature_point_curve_on_surface_geometry.h:55
QuadraturePointCurveOnSurfaceGeometry(QuadraturePointCurveOnSurfaceGeometry const &rOther)
Copy constructor.
Definition: quadrature_point_curve_on_surface_geometry.h:109
A sinlge quadrature point, that can be used for geometries without a predefined integration scheme,...
Definition: quadrature_point_geometry.h:44
GeometryType & GetGeometryParent(IndexType Index) const override
Definition: quadrature_point_geometry.h:310
const Matrix & ShapeFunctionsValues() const
Definition: geometry.h:3393
Vector & DeterminantOfJacobian(Vector &rResult) const
Definition: geometry.h:3154
QuadraturePointGeometry & operator=(const QuadraturePointGeometry &rOther)
Assignment operator.
Definition: quadrature_point_geometry.h:220
JacobiansType & InverseOfJacobian(JacobiansType &rResult) const
Definition: geometry.h:3262
const ShapeFunctionsGradientsType & ShapeFunctionsLocalGradients() const
Definition: geometry.h:3539
JacobiansType & Jacobian(JacobiansType &rResult) const
using base class functions
Definition: geometry.h:2901
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_DEBUG_ERROR_IF(conditional)
Definition: exception.h:171
REF: G. R. Cowper, GAUSSIAN QUADRATURE FORMULAS FOR TRIANGLES.
Definition: mesh_condition.cpp:21
TExpressionType::data_type norm_2(AMatrix::MatrixExpression< TExpressionType, TCategory > const &TheExpression)
Definition: amatrix_interface.h:625
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