12#ifndef DUMUX_DISCRETIZATION_HYBRID_CVFE_FLUXVARIABLES_CACHE_HH
13#define DUMUX_DISCRETIZATION_HYBRID_CVFE_FLUXVARIABLES_CACHE_HH
17#include <dune/common/fvector.hh>
18#include <dumux/common/typetraits/localdofs_.hh>
19#include <dumux/common/concepts/ipdata_.hh>
30template<
class Scalar,
class Gr
idGeometry >
33 using GridView =
typename GridGeometry::GridView;
34 using Element =
typename GridView::template Codim<0>::Entity;
35 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
36 using LocalPosition =
typename Element::Geometry::LocalCoordinate;
38 using FVElementGeometry =
typename GridGeometry::LocalView;
39 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
41 static const int dim = GridView::dimension;
42 static const int dimWorld = GridView::dimensionworld;
44 using CoordScalar =
typename GridView::ctype;
45 using FeLocalBasis =
typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
46 using ShapeJacobian =
typename FeLocalBasis::Traits::JacobianType;
47 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
48 using JacobianInverseTransposed =
typename Element::Geometry::JacobianInverseTransposed;
57 template<
class Problem,
class ElementVolumeVariables >
59 const Element& element,
60 const FVElementGeometry& fvGeometry,
61 const ElementVolumeVariables& elemVolVars,
62 const GlobalPosition& globalPos)
64 update_(fvGeometry, fvGeometry.elementGeometry().local(globalPos), globalPos);
68 template<
class Problem,
class ElementVolumeVariables, Concept::IpData IpData >
70 const Element& element,
71 const FVElementGeometry& fvGeometry,
72 const ElementVolumeVariables& elemVolVars,
75 update_(fvGeometry, ipData.local(), ipData.global());
79 const GlobalPosition&
ipGlobal()
const {
return ipGlobal_; }
81 const LocalPosition&
ipLocal()
const {
return ipLocal_; }
83 const std::vector<ShapeJacobian>&
shapeJacobian()
const {
return shapeJacobian_; }
85 const std::vector<ShapeValue>&
shapeValues()
const {
return shapeValues_; }
87 const JacobianInverseTransposed&
jacInvT()
const {
return jacInvT_; }
89 const GlobalPosition&
gradN(
unsigned int scvIdxInElement)
const {
return gradN_[scvIdxInElement]; }
93 void update_(
const FVElementGeometry& fvGeometry,
94 const LocalPosition& localPos,
95 const GlobalPosition& globalPos)
97 const auto& geometry = fvGeometry.elementGeometry();
98 const auto& localBasis = fvGeometry.feLocalBasis();
101 ipGlobal_ = globalPos;
103 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal_);
104 localBasis.evaluateJacobian(ipLocal_, shapeJacobian_);
105 localBasis.evaluateFunction(ipLocal_, shapeValues_);
108 gradN_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry));
109 for (
const auto& localDof:
localDofs(fvGeometry))
110 jacInvT_.mv(shapeJacobian_[localDof.index()][0], gradN_[localDof.index()]);
113 std::vector<GlobalPosition> gradN_;
114 std::vector<ShapeJacobian> shapeJacobian_;
115 std::vector<ShapeValue> shapeValues_;
116 JacobianInverseTransposed jacInvT_;
117 GlobalPosition ipGlobal_;
118 LocalPosition ipLocal_;
Flux variables cache class for control-volume finite element schemes. For control-volume finite eleme...
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:32
static bool constexpr isSolDependent
whether the cache needs an update when the solution changes
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:54
const LocalPosition & ipLocal() const
returns the local position for which this cache has been updated
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:81
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
returns the shape function gradients in global coordinates at the interpolation point
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:89
const std::vector< ShapeJacobian > & shapeJacobian() const
returns the shape function gradients in local coordinates at the interpolation point
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:83
const JacobianInverseTransposed & jacInvT() const
returns inverse transposed jacobian at the interpolation point
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:87
const std::vector< ShapeValue > & shapeValues() const
returns the shape function values at the interpolation point
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:85
const GlobalPosition & ipGlobal() const
returns the global position for which this cache has been updated
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:79
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/hybrid/fluxvariablescache.hh:58
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const IpData &ipData)
update the cache for interpolation point data
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:69
typename GridGeometry::ScvfQuadratureRule ScvfQuadratureRule
Definition: discretization/cvfe/hybrid/fluxvariablescache.hh:51
Class representing dofs on elements for control-volume finite element schemes.
auto localDofs(const FVElementGeometry &fvGeometry)
range over local dofs
Definition: localdof.hh:50