15 #if !defined(KRATOS_QUADRATURE_POINT_CURVE_ON_SURFACE_GEOMETRY_H_INCLUDED )
16 #define KRATOS_QUADRATURE_POINT_CURVE_ON_SURFACE_GEOMETRY_H_INCLUDED
40 template<
class TPo
intType>
84 double LocalTangentsU,
85 double LocalTangentsV)
86 :
BaseType(ThisPoints, ThisGeometryShapeFunctionContainer)
87 , mLocalTangentsU(LocalTangentsU)
88 , mLocalTangentsV(LocalTangentsV)
96 double LocalTangentsU,
97 double LocalTangentsV,
99 :
BaseType(ThisPoints, ThisGeometryShapeFunctionContainer, pGeometryParent)
100 , mLocalTangentsU(LocalTangentsU)
101 , mLocalTangentsV(LocalTangentsV)
111 , mLocalTangentsU(rOther.mLocalTangentsU)
112 , mLocalTangentsV(rOther.mLocalTangentsV)
125 mLocalTangentsU = rOther.mLocalTangentsU;
126 mLocalTangentsV = rOther.mLocalTangentsV;
140 if (rVariable == LOCAL_TANGENT)
142 mLocalTangentsU = rInput[0];
143 mLocalTangentsV = rInput[1];
152 if (rVariable == LOCAL_TANGENT)
154 rOutput[0] = mLocalTangentsU;
155 rOutput[1] = mLocalTangentsV;
163 Vector& rOutput)
const override
165 if (rVariable == DETERMINANTS_OF_JACOBIAN_PARENT)
186 <<
"Trying to access Normal of QuadraturePointCurveOnSurface "
187 <<
"with an integration point index != 0." << std::endl;
190 this->
Jacobian(J, IntegrationPointIndex, ThisMethod);
214 <<
"Trying to access DeterminantOfJacobian of QuadraturePointCurveOnSurface "
215 <<
"with an integration point index != 0." << std::endl;
218 this->
Jacobian(J, IntegrationPointIndex, ThisMethod);
223 return norm_2(a_1 * mLocalTangentsU + a_2 * mLocalTangentsV);
234 if (rResult.size() != 1)
251 if (rResult.size() != 1)
281 std::string
Info()
const override
283 return "Quadrature point for a curve on surface.";
289 rOStream <<
"Quadrature point for a curve on surface.";
302 double mLocalTangentsU;
303 double mLocalTangentsV;
317 void save(
Serializer& rSerializer )
const override
320 rSerializer.
save(
"LocalTangentsU", mLocalTangentsU);
321 rSerializer.
save(
"LocalTangentsV", mLocalTangentsV);
327 rSerializer.
load(
"LocalTangentsU", mLocalTangentsU);
328 rSerializer.
load(
"LocalTangentsV", mLocalTangentsV);
KratosGeometryType
Definition: geometry_data.h:110
@ Kratos_Quadrature_Point_Curve_On_Surface_Geometry
IntegrationMethod
Definition: geometry_data.h:76
KratosGeometryFamily
Definition: geometry_data.h:91
@ Kratos_Quadrature_Geometry
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