26#ifndef DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
27#define DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
29#include <dune/localfunctions/lagrange/pqkfactory.hh>
43template<
class TypeTag>
52 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
53 using Element =
typename GridView::template Codim<0>::Entity;
54 using GridIndexType =
typename GridView::IndexSet::IndexType;
55 using Stencil = std::vector<GridIndexType>;
56 using TransmissibilityVector = std::vector<GridIndexType>;
58 using CoordScalar =
typename GridView::ctype;
59 static const int dim = GridView::dimension;
60 static const int dimWorld = GridView::dimensionworld;
62 using FeCache = Dune::PQkLocalFiniteElementCache<CoordScalar, Scalar, dim, 1>;
63 using FeLocalBasis =
typename FeCache::FiniteElementType::Traits::LocalBasisType;
64 using ShapeJacobian =
typename FeLocalBasis::Traits::JacobianType;
65 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
66 using JacobianInverseTransposed =
typename Element::Geometry::JacobianInverseTransposed;
67 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
72 const Element& element,
73 const FVElementGeometry& fvGeometry,
74 const ElementVolumeVariables& elemVolVars,
75 const SubControlVolumeFace& scvf)
77 const auto geometry = element.geometry();
78 const auto& localBasis = fvGeometry.feLocalBasis();
81 std::vector<ShapeValue> shapeVals;
82 const auto ipLocal = geometry.local(scvf.ipGlobal());
83 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal);
84 localBasis.evaluateJacobian(ipLocal, shapeJacobian_);
85 localBasis.evaluateFunction(ipLocal, shapeVals);
88 shapeValues_.resize(fvGeometry.numScv(), 0.0);
89 if (!scvf.isOnFracture())
90 std::copy(shapeVals.begin(), shapeVals.end(), shapeValues_.begin());
93 const auto thisFacetIdx = scvf.facetIndexInElement();
94 for (
const auto& scv: scvs(fvGeometry))
95 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
96 shapeValues_[scv.indexInElement()] = shapeVals[scv.localDofIndex()];
100 gradN_.resize(fvGeometry.numScv(), GlobalPosition(0.0));
101 if (!scvf.isOnFracture())
103 for (
const auto& scv: scvs(fvGeometry))
104 if (!scv.isOnFracture())
105 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
109 const auto thisFacetIdx = scvf.facetIndexInElement();
112 std::vector<unsigned int> facetLocalDofs;
113 for (
const auto& scv : scvs(fvGeometry))
114 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
115 facetLocalDofs.push_back(scv.localDofIndex());
117 for (
const auto& scv: scvs(fvGeometry))
120 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
121 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
124 else if (!scv.isOnFracture()
125 && std::find( facetLocalDofs.begin(),
126 facetLocalDofs.end(),
127 scv.localDofIndex() ) == facetLocalDofs.end())
129 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
136 const std::vector<ShapeJacobian>&
shapeJacobian()
const {
return shapeJacobian_; }
138 const std::vector<ShapeValue>&
shapeValues()
const {
return shapeValues_; }
140 const JacobianInverseTransposed&
jacInvT()
const {
return jacInvT_; }
142 const GlobalPosition&
gradN(
unsigned int scvIdxInElement)
const {
return gradN_[scvIdxInElement]; }
145 std::vector<GlobalPosition> gradN_;
146 std::vector<ShapeJacobian> shapeJacobian_;
147 std::vector<ShapeValue> shapeValues_;
148 JacobianInverseTransposed jacInvT_;
The available discretization methods in Dumux.
Classes related to flux variables caching.
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
We only store discretization-related quantities for the box method. However, we cannot reuse the cach...
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:45
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
Returns the shape value gradients corresponding to an scv.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:142
const std::vector< ShapeValue > & shapeValues() const
Returns the shape values for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:138
const std::vector< ShapeJacobian > & shapeJacobian() const
Returns the Jacobian of the shape functions at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:136
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:71
const JacobianInverseTransposed & jacInvT() const
Returns the shape value gradients for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:140
Declares all properties used in Dumux.