12#ifndef DUMUX_DISCRETIZATION_CVFE_FLUXVARIABLES_CACHE_HH
13#define DUMUX_DISCRETIZATION_CVFE_FLUXVARIABLES_CACHE_HH
15#include <dune/common/fvector.hh>
16#include <dumux/common/typetraits/localdofs_.hh>
28template<
class Scalar,
class Gr
idGeometry >
31 using GridView =
typename GridGeometry::GridView;
32 using Element =
typename GridView::template Codim<0>::Entity;
33 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
34 using LocalPosition =
typename Element::Geometry::LocalCoordinate;
36 using FVElementGeometry =
typename GridGeometry::LocalView;
37 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
39 static const int dim = GridView::dimension;
40 static const int dimWorld = GridView::dimensionworld;
42 using CoordScalar =
typename GridView::ctype;
43 using FeLocalBasis =
typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
44 using ShapeJacobian =
typename FeLocalBasis::Traits::JacobianType;
45 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
46 using JacobianInverseTransposed =
typename Element::Geometry::JacobianInverseTransposed;
55 template<
class Problem,
class ElementVolumeVariables >
57 const Element& element,
58 const FVElementGeometry& fvGeometry,
59 const ElementVolumeVariables& elemVolVars,
60 const SubControlVolumeFace& scvf)
62 update(problem, element, fvGeometry, elemVolVars, scvf.ipGlobal());
66 template<
class Problem,
class ElementVolumeVariables >
68 const Element& element,
69 const FVElementGeometry& fvGeometry,
70 const ElementVolumeVariables& elemVolVars,
71 const GlobalPosition& globalPos)
73 const auto geometry = element.geometry();
74 const auto& localBasis = fvGeometry.feLocalBasis();
77 ipGlobal_ = globalPos;
78 ipLocal_ = geometry.local(globalPos);
79 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal_);
80 localBasis.evaluateJacobian(ipLocal_, shapeJacobian_);
81 localBasis.evaluateFunction(ipLocal_, shapeValues_);
84 gradN_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
85 for (
const auto& localDof:
localDofs(fvGeometry))
86 jacInvT_.mv(shapeJacobian_[localDof.index()][0], gradN_[localDof.index()]);
90 const GlobalPosition&
ipGlobal()
const {
return ipGlobal_; }
92 const GlobalPosition&
ipLocal()
const {
return ipLocal_; }
96 const std::vector<ShapeJacobian>&
shapeJacobian()
const {
return shapeJacobian_; }
98 const std::vector<ShapeValue>&
shapeValues()
const {
return shapeValues_; }
100 const JacobianInverseTransposed&
jacInvT()
const {
return jacInvT_; }
102 const GlobalPosition&
gradN(
unsigned int scvIdxInElement)
const {
return gradN_[scvIdxInElement]; }
105 std::vector<GlobalPosition> gradN_;
106 std::vector<ShapeJacobian> shapeJacobian_;
107 std::vector<ShapeValue> shapeValues_;
108 JacobianInverseTransposed jacInvT_;
109 GlobalPosition ipGlobal_;
110 LocalPosition ipLocal_;
An interpolation point related to an element that includes global and local positions.
Definition: cvfe/interpolationpointdata.hh:25
Flux variables cache class for control-volume finite element schemes. For control-volume finite eleme...
Definition: discretization/cvfe/fluxvariablescache.hh:30
const JacobianInverseTransposed & jacInvT() const
returns inverse transposed jacobian at the interpolation point
Definition: discretization/cvfe/fluxvariablescache.hh:100
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
update the cache for an scvf
Definition: discretization/cvfe/fluxvariablescache.hh:56
const std::vector< ShapeValue > & shapeValues() const
returns the shape function values at the interpolation point
Definition: discretization/cvfe/fluxvariablescache.hh:98
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
returns the shape function gradients in global coordinates at the interpolation point
Definition: discretization/cvfe/fluxvariablescache.hh:102
const GlobalPosition & ipGlobal() const
returns the global position for which this cache has been updated
Definition: discretization/cvfe/fluxvariablescache.hh:90
const GlobalPosition & ipLocal() const
returns the local position for which this cache has been updated
Definition: discretization/cvfe/fluxvariablescache.hh:92
IpData ipData() const
returns the ipData for which this cache has been updated
Definition: discretization/cvfe/fluxvariablescache.hh:94
static bool constexpr isSolDependent
whether the cache needs an update when the solution changes
Definition: discretization/cvfe/fluxvariablescache.hh:52
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const GlobalPosition &globalPos)
update the cache for a given global position
Definition: discretization/cvfe/fluxvariablescache.hh:67
const std::vector< ShapeJacobian > & shapeJacobian() const
returns the shape function gradients in local coordinates at the interpolation point
Definition: discretization/cvfe/fluxvariablescache.hh:96
Classes representing interpolation point data for control-volume finite element schemes.
Class representing dofs on elements for control-volume finite element schemes.
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50