14#ifndef DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
15#define DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
17#include <dune/common/fvector.hh>
31template<
class TypeTag>
37 using GridView =
typename GridGeometry::GridView;
39 using FVElementGeometry =
typename GridGeometry::LocalView;
41 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
42 using Element =
typename GridView::template Codim<0>::Entity;
43 using GridIndexType =
typename GridView::IndexSet::IndexType;
44 using Stencil = std::vector<GridIndexType>;
45 using TransmissibilityVector = std::vector<GridIndexType>;
47 using CoordScalar =
typename GridView::ctype;
48 static const int dim = GridView::dimension;
49 static const int dimWorld = GridView::dimensionworld;
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;
55 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
61 const Element& element,
62 const FVElementGeometry& fvGeometry,
63 const ElementVolumeVariables& elemVolVars,
64 const SubControlVolumeFace& scvf)
66 const auto geometry = element.geometry();
67 const auto& localBasis = fvGeometry.feLocalBasis();
70 std::vector<ShapeValue> shapeVals;
71 const auto ipLocal = geometry.local(scvf.ipGlobal());
72 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal);
73 localBasis.evaluateJacobian(ipLocal, shapeJacobian_);
74 localBasis.evaluateFunction(ipLocal, shapeVals);
77 shapeValues_.resize(fvGeometry.numScv(), 0.0);
78 if (!scvf.isOnFracture())
79 std::copy(shapeVals.begin(), shapeVals.end(), shapeValues_.begin());
82 const auto thisFacetIdx = scvf.facetIndexInElement();
83 for (
const auto& scv: scvs(fvGeometry))
84 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
85 shapeValues_[scv.indexInElement()] = shapeVals[scv.localDofIndex()];
89 gradN_.resize(fvGeometry.numScv(), GlobalPosition(0.0));
90 if (!scvf.isOnFracture())
92 for (
const auto& scv: scvs(fvGeometry))
93 if (!scv.isOnFracture())
94 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
98 const auto thisFacetIdx = scvf.facetIndexInElement();
101 std::vector<unsigned int> facetLocalDofs;
102 for (
const auto& scv : scvs(fvGeometry))
103 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
104 facetLocalDofs.push_back(scv.localDofIndex());
106 for (
const auto& scv: scvs(fvGeometry))
109 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
110 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
113 else if (!scv.isOnFracture()
114 && std::find( facetLocalDofs.begin(),
115 facetLocalDofs.end(),
116 scv.localDofIndex() ) == facetLocalDofs.end())
118 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
125 const std::vector<ShapeJacobian>&
shapeJacobian()
const {
return shapeJacobian_; }
127 const std::vector<ShapeValue>&
shapeValues()
const {
return shapeValues_; }
129 const JacobianInverseTransposed&
jacInvT()
const {
return jacInvT_; }
131 const GlobalPosition&
gradN(
unsigned int scvIdxInElement)
const {
return gradN_[scvIdxInElement]; }
134 std::vector<GlobalPosition> gradN_;
135 std::vector<ShapeJacobian> shapeJacobian_;
136 std::vector<ShapeValue> shapeValues_;
137 JacobianInverseTransposed jacInvT_;
We only store discretization-related quantities for the box method. However, we cannot reuse the cach...
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:33
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
Returns the shape value gradients corresponding to an scv.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:131
static constexpr bool isSolDependent
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:58
const std::vector< ShapeValue > & shapeValues() const
Returns the shape values for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:127
const std::vector< ShapeJacobian > & shapeJacobian() const
Returns the Jacobian of the shape functions at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:125
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:60
const JacobianInverseTransposed & jacInvT() const
Returns the shape value gradients for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:129
Defines all properties used in Dumux.
Classes related to flux variables caching.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.