12#ifndef DUMUX_DISCRETIZATION_FE_INTEGRATION_POINT_DATA_HH
13#define DUMUX_DISCRETIZATION_FE_INTEGRATION_POINT_DATA_HH
19template<
class GlobalPosition,
class LocalBasis>
22 using LocalPosition =
typename LocalBasis::Traits::DomainType;
23 using RangeType =
typename LocalBasis::Traits::RangeType;
24 using JacobianType =
typename LocalBasis::Traits::JacobianType;
26 using ShapeValues = std::vector<RangeType>;
27 using ShapeGradients = std::vector<GlobalPosition>;
34 template<
class Geometry>
36 const LocalPosition&
local,
37 const LocalBasis& localBasis)
41 auto numLocalDofs = localBasis.size();
44 shapeValues_.resize(numLocalDofs);
45 localBasis.evaluateFunction(
local, shapeValues_);
48 std::vector<JacobianType> shapeGrads(numLocalDofs);
49 localBasis.evaluateJacobian(
local, shapeGrads);
52 const auto jacInvT = geometry.jacobianInverseTransposed(
local);
53 shapeGradients_.resize(numLocalDofs, GlobalPosition(0.0));
54 for (
unsigned int i = 0; i < numLocalDofs; ++i)
55 jacInvT.umv(shapeGrads[i][0], shapeGradients_[i]);
60 {
return shapeValues_; }
64 {
return shapeValues_[i]; }
68 {
return shapeGradients_; }
71 const GlobalPosition&
gradN(
int i)
const
72 {
return shapeGradients_[i]; }
75 const LocalPosition&
local()
const
84 GlobalPosition global_;
86 ShapeValues shapeValues_;
87 ShapeGradients shapeGradients_;
94template<
class GlobalPosition,
class LocalBasis,
class BoundaryFlag>
98 using LocalPosition =
typename LocalBasis::Traits::DomainType;
104 template<
class Geometry>
106 const LocalPosition&
local,
107 const LocalBasis& localBasis,
108 const GlobalPosition& n,
110 :
ParentType(geometry,
local, localBasis), normal_(n), boundaryFlag_(bFlag)
119 {
return boundaryFlag_.
get(); }
122 const GlobalPosition& normal_;
Boundary flag to store e.g. in sub control volume faces.
Definition: boundaryflag.hh:58
std::size_t value_type
Definition: boundaryflag.hh:28
value_type get() const
Definition: boundaryflag.hh:41
Interpolation point data related to a face of an element.
Definition: fem/interpolationpointdata.hh:96
BoundaryFlag::value_type boundaryFlag() const
Return the boundary flag.
Definition: fem/interpolationpointdata.hh:118
FEFaceInterpolationPointData()=delete
FEFaceInterpolationPointData(const Geometry &geometry, const LocalPosition &local, const LocalBasis &localBasis, const GlobalPosition &n, const BoundaryFlag &bFlag)
Definition: fem/interpolationpointdata.hh:105
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector at the quadrature point.
Definition: fem/interpolationpointdata.hh:114
Definition: fem/interpolationpointdata.hh:21
const ShapeValues & shapeValues() const
The shape values at the quadrature point.
Definition: fem/interpolationpointdata.hh:59
const LocalPosition & local() const
The local position of the quadrature point.
Definition: fem/interpolationpointdata.hh:75
const GlobalPosition & global() const
The global position of the quadrature point.
Definition: fem/interpolationpointdata.hh:79
FEInterpolationPointData()=delete
const RangeType & shapeValue(int i) const
The shape value of a local dof at the quadrature point.
Definition: fem/interpolationpointdata.hh:63
FEInterpolationPointData(const Geometry &geometry, const LocalPosition &local, const LocalBasis &localBasis)
Definition: fem/interpolationpointdata.hh:35
const ShapeGradients & shapeGradients() const
The shape value gradients at the quadrature point.
Definition: fem/interpolationpointdata.hh:67
const GlobalPosition & gradN(int i) const
The shape value gradient of a local dof at the quadrature point.
Definition: fem/interpolationpointdata.hh:71