12#ifndef DUMUX_DISCRETIZATION_CVFE_FLUXVARIABLES_CACHE_HH
13#define DUMUX_DISCRETIZATION_CVFE_FLUXVARIABLES_CACHE_HH
15#include <dune/common/fvector.hh>
25template<
class Scalar,
class Gr
idGeometry >
28 using GridView =
typename GridGeometry::GridView;
29 using Element =
typename GridView::template Codim<0>::Entity;
30 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
32 using FVElementGeometry =
typename GridGeometry::LocalView;
33 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
35 static const int dim = GridView::dimension;
36 static const int dimWorld = GridView::dimensionworld;
38 using CoordScalar =
typename GridView::ctype;
39 using FeLocalBasis =
typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
40 using ShapeJacobian =
typename FeLocalBasis::Traits::JacobianType;
41 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
42 using JacobianInverseTransposed =
typename Element::Geometry::JacobianInverseTransposed;
49 template<
class Problem,
class ElementVolumeVariables >
51 const Element& element,
52 const FVElementGeometry& fvGeometry,
53 const ElementVolumeVariables& elemVolVars,
54 const SubControlVolumeFace& scvf)
56 update(problem, element, fvGeometry, elemVolVars, scvf.ipGlobal());
60 template<
class Problem,
class ElementVolumeVariables >
62 const Element& element,
63 const FVElementGeometry& fvGeometry,
64 const ElementVolumeVariables& elemVolVars,
65 const GlobalPosition& globalPos)
67 const auto geometry = element.geometry();
68 const auto& localBasis = fvGeometry.feLocalBasis();
71 ipGlobal_ = globalPos;
72 const auto ipLocal = geometry.local(ipGlobal_);
73 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal);
74 localBasis.evaluateJacobian(ipLocal, shapeJacobian_);
75 localBasis.evaluateFunction(ipLocal, shapeValues_);
78 gradN_.resize(fvGeometry.numScv());
79 for (
const auto& scv: scvs(fvGeometry))
80 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
84 const GlobalPosition&
ipGlobal()
const {
return ipGlobal_; }
86 const std::vector<ShapeJacobian>&
shapeJacobian()
const {
return shapeJacobian_; }
88 const std::vector<ShapeValue>&
shapeValues()
const {
return shapeValues_; }
90 const JacobianInverseTransposed&
jacInvT()
const {
return jacInvT_; }
92 const GlobalPosition&
gradN(
unsigned int scvIdxInElement)
const {
return gradN_[scvIdxInElement]; }
95 GlobalPosition ipGlobal_;
96 std::vector<GlobalPosition> gradN_;
97 std::vector<ShapeJacobian> shapeJacobian_;
98 std::vector<ShapeValue> shapeValues_;
99 JacobianInverseTransposed jacInvT_;
Flux variables cache class for control-volume finite element schemes. For control-volume finite eleme...
Definition: discretization/cvfe/fluxvariablescache.hh:27
const JacobianInverseTransposed & jacInvT() const
returns inverse transposed jacobian at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:90
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:50
const std::vector< ShapeValue > & shapeValues() const
returns the shape function values at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:88
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
returns the shape function gradients in global coordinates at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:92
const GlobalPosition & ipGlobal() const
returns the global position for which this cache has been updated
Definition: discretization/cvfe/fluxvariablescache.hh:84
static bool constexpr isSolDependent
whether the cache needs an update when the solution changes
Definition: discretization/cvfe/fluxvariablescache.hh:46
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:61
const std::vector< ShapeJacobian > & shapeJacobian() const
returns the shape function gradients in local coordinates at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:86