24#ifndef DUMUX_DISCRETIZATION_CVFE_FLUXVARIABLES_CACHE_HH
25#define DUMUX_DISCRETIZATION_CVFE_FLUXVARIABLES_CACHE_HH
27#include <dune/common/fvector.hh>
37template<
class Scalar,
class Gr
idGeometry >
40 using GridView =
typename GridGeometry::GridView;
41 using Element =
typename GridView::template Codim<0>::Entity;
42 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
44 using FVElementGeometry =
typename GridGeometry::LocalView;
45 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
47 static const int dim = GridView::dimension;
48 static const int dimWorld = GridView::dimensionworld;
50 using CoordScalar =
typename GridView::ctype;
51 using FeLocalBasis =
typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
52 using ShapeJacobian =
typename FeLocalBasis::Traits::JacobianType;
53 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
54 using JacobianInverseTransposed =
typename Element::Geometry::JacobianInverseTransposed;
59 template<
class Problem,
class ElementVolumeVariables >
61 const Element& element,
62 const FVElementGeometry& fvGeometry,
63 const ElementVolumeVariables& elemVolVars,
64 const SubControlVolumeFace& scvf)
66 update(problem, element, fvGeometry, elemVolVars, scvf.ipGlobal());
70 template<
class Problem,
class ElementVolumeVariables >
72 const Element& element,
73 const FVElementGeometry& fvGeometry,
74 const ElementVolumeVariables& elemVolVars,
75 const GlobalPosition& globalPos)
77 const auto geometry = element.geometry();
78 const auto& localBasis = fvGeometry.feLocalBasis();
81 ipGlobal_ = globalPos;
82 const auto ipLocal = geometry.local(ipGlobal_);
83 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal);
84 localBasis.evaluateJacobian(ipLocal, shapeJacobian_);
85 localBasis.evaluateFunction(ipLocal, shapeValues_);
88 gradN_.resize(fvGeometry.numScv());
89 for (
const auto& scv: scvs(fvGeometry))
90 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
94 const GlobalPosition&
ipGlobal()
const {
return ipGlobal_; }
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 GlobalPosition ipGlobal_;
106 std::vector<GlobalPosition> gradN_;
107 std::vector<ShapeJacobian> shapeJacobian_;
108 std::vector<ShapeValue> shapeValues_;
109 JacobianInverseTransposed jacInvT_;
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
Flux variables cache class for control-volume finite element schemes. For control-volume finite eleme...
Definition: discretization/cvfe/fluxvariablescache.hh:39
const JacobianInverseTransposed & jacInvT() const
returns inverse transposed jacobian at the integration 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:60
const std::vector< ShapeValue > & shapeValues() const
returns the shape function values at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:98
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
returns the shape function gradients in global coordinates at the integration 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:94
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:71
const std::vector< ShapeJacobian > & shapeJacobian() const
returns the shape function gradients in local coordinates at the integration point
Definition: discretization/cvfe/fluxvariablescache.hh:96