12#ifndef DUMUX_CVFE_IP_DATA_HH
13#define DUMUX_CVFE_IP_DATA_HH
19#include <dumux/common/typetraits/localdofs_.hh>
21#include <dumux/common/concepts/ipdata_.hh>
29template<
class LocalPos,
class GlobalPos>
37 : local_(std::move(localPos)), global_(std::move(pos)) {}
39 : local_(localPos), global_(pos) {}
58template<
class LocalPosition,
class GlobalPosition,
class LocalIndex>
64 :
ParentType(localPos, pos), localDofIndex_(index) {}
66 :
ParentType(localPos, pos), localDofIndex_(index) {}
70 {
return localDofIndex_; }
73 LocalIndex localDofIndex_;
80template<
class LocalMapping,
class GlobalPos>
84 using LocalPosition = std::invoke_result_t<LocalMapping, const GlobalPos&>;
96 {
return localMapping_(global_); }
99 LocalMapping localMapping_;
107template<
class BaseClass,
class LocalIndex>
111 using GlobalPosition = std::remove_cvref_t<decltype(std::declval<BaseClass>().global())>;
112 using LocalPosition = std::remove_cvref_t<decltype(std::declval<BaseClass>().local())>;
114 template<
class... Args>
116 : BaseClass(std::forward<Args>(args)...), normal_(std::move(n)), scvfIndex_(index) {}
118 template<
class... Args>
120 : BaseClass(std::forward<Args>(args)...), normal_(n), scvfIndex_(index) {}
128 {
return scvfIndex_; }
132 LocalIndex scvfIndex_;
139template<
class IpData>
143 template<
class... Args>
145 : IpData(std::forward<Args>(args)...), qpIndex_(qpIdx) {}
152 std::size_t qpIndex_;
159template<
class Gr
idGeometry >
162 using GridView =
typename GridGeometry::GridView;
163 using Element =
typename GridView::template Codim<0>::Entity;
164 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
166 using LocalBasis =
typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
167 using RangeType =
typename LocalBasis::Traits::RangeType;
168 using JacobianType =
typename LocalBasis::Traits::JacobianType;
170 using Gradients = std::vector<GlobalPosition>;
172 using FVElementGeometry =
typename GridGeometry::LocalView;
179 template<
class Problem,
class ElementVariables >
181 const Element& element,
182 const FVElementGeometry& fvGeometry,
183 const ElementVariables& elemVars,
184 const GlobalPosition& globalPos)
186 update_(fvGeometry, fvGeometry.elementGeometry().local(globalPos));
190 template<
class Problem,
class ElementVariables, Concept::IpData IpData >
192 const Element& element,
193 const FVElementGeometry& fvGeometry,
194 const ElementVariables& elemVars,
195 const IpData& ipData)
197 update_(fvGeometry, ipData.local());
201 const std::vector<JacobianType>&
shapeJacobian()
const {
return shapeJacobian_; }
203 const std::vector<RangeType>&
shapeValues()
const {
return shapeValues_; }
205 const GlobalPosition&
gradN(
unsigned int localIdx)
const {
return gradN_[localIdx]; }
209 template<
class LocalPosition>
210 void update_(
const FVElementGeometry& fvGeometry,
211 const LocalPosition& localPos)
213 const auto& geometry = fvGeometry.elementGeometry();
214 const auto& localBasis = fvGeometry.feLocalBasis();
217 const auto jacInvT = geometry.jacobianInverseTransposed(localPos);
218 localBasis.evaluateJacobian(localPos, shapeJacobian_);
219 localBasis.evaluateFunction(localPos, shapeValues_);
222 gradN_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
223 for (
const auto& localDof:
localDofs(fvGeometry))
224 jacInvT.mv(shapeJacobian_[localDof.index()][0], gradN_[localDof.index()]);
228 std::vector<JacobianType> shapeJacobian_;
229 std::vector<RangeType> shapeValues_;
An interpolation point related to a face of an element.
Definition: cvfe/interpolationpointdata.hh:109
LocalIndex scvfIndex() const
The local index of an scvf.
Definition: cvfe/interpolationpointdata.hh:127
FaceInterpolationPointData(const GlobalPosition &n, LocalIndex index, Args &&... args)
Definition: cvfe/interpolationpointdata.hh:119
FaceInterpolationPointData(GlobalPosition &&n, LocalIndex index, Args &&... args)
Definition: cvfe/interpolationpointdata.hh:115
const GlobalPosition & unitOuterNormal() const
The unit outer normal vector at the quadrature point.
Definition: cvfe/interpolationpointdata.hh:123
std::remove_cvref_t< decltype(std::declval< BaseClass >().local())> LocalPosition
Definition: cvfe/interpolationpointdata.hh:112
std::remove_cvref_t< decltype(std::declval< BaseClass >().global())> GlobalPosition
Definition: cvfe/interpolationpointdata.hh:111
Wraps interpolation point data and adds a quadrature point index for use in quadrature loops.
Definition: cvfe/interpolationpointdata.hh:141
std::size_t qpIndex() const
The quadrature point index.
Definition: cvfe/interpolationpointdata.hh:148
IndexedQuadratureInterpolationPointData(std::size_t qpIdx, Args &&... args)
Definition: cvfe/interpolationpointdata.hh:144
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:31
const LocalPosition & local() const
The local position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:46
GlobalPos GlobalPosition
Definition: cvfe/interpolationpointdata.hh:34
const GlobalPosition & global() const
The global position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:42
InterpolationPointData(LocalPosition &&localPos, GlobalPosition &&pos)
Definition: cvfe/interpolationpointdata.hh:36
LocalPos LocalPosition
Definition: cvfe/interpolationpointdata.hh:33
InterpolationPointData(const LocalPosition &localPos, const GlobalPosition &pos)
Definition: cvfe/interpolationpointdata.hh:38
An interpolation point related to a global position of an element, giving its local positions by a ma...
Definition: cvfe/interpolationpointdata.hh:82
InterpolationPointDataLocalMapping(LocalMapping &&mapping, GlobalPosition &&pos)
Definition: cvfe/interpolationpointdata.hh:87
const LocalPosition local() const
The local position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:95
InterpolationPointDataLocalMapping(LocalMapping &&mapping, const GlobalPosition &pos)
Definition: cvfe/interpolationpointdata.hh:88
std::invoke_result_t< LocalMapping, const GlobalPos & > LocalPosition
Definition: cvfe/interpolationpointdata.hh:84
const GlobalPosition & global() const
The global position of the quadrature point.
Definition: cvfe/interpolationpointdata.hh:91
GlobalPos GlobalPosition
Definition: cvfe/interpolationpointdata.hh:85
Interpolation point data related to a local basis.
Definition: cvfe/interpolationpointdata.hh:161
static constexpr bool isSolDependent
whether the cache needs an update when the solution changes
Definition: cvfe/interpolationpointdata.hh:176
const GlobalPosition & gradN(unsigned int localIdx) const
returns the shape function gradients in global coordinates at the interpolation point
Definition: cvfe/interpolationpointdata.hh:205
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const GlobalPosition &globalPos)
update the cache for a given global position
Definition: cvfe/interpolationpointdata.hh:180
const std::vector< JacobianType > & shapeJacobian() const
returns the shape function gradients in local coordinates at the interpolation point
Definition: cvfe/interpolationpointdata.hh:201
const std::vector< RangeType > & shapeValues() const
returns the shape function values at the interpolation point
Definition: cvfe/interpolationpointdata.hh:203
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVariables &elemVars, const IpData &ipData)
update the cache for interpolation point data
Definition: cvfe/interpolationpointdata.hh:191
An interpolation point related to a localDof of an element, giving its global and local positions.
Definition: cvfe/interpolationpointdata.hh:60
LocalDofInterpolationPointData(const LocalPosition &localPos, const GlobalPosition &pos, LocalIndex index)
Definition: cvfe/interpolationpointdata.hh:65
LocalDofInterpolationPointData(LocalPosition &&localPos, GlobalPosition &&pos, LocalIndex index)
Definition: cvfe/interpolationpointdata.hh:63
LocalIndex localDofIndex() const
The local index of the corresponding dof.
Definition: cvfe/interpolationpointdata.hh:69
Class representing dofs on elements for control-volume finite element schemes.
Definition: cvfe/interpolationpointdata.hh:23
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50