26#ifndef DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
27#define DUMUX_POROUSMEDIUM_BOXDFM_FLUXVARIABLESCACHE_HH
41template<
class TypeTag>
47 using GridView =
typename GridGeometry::GridView;
49 using FVElementGeometry =
typename GridGeometry::LocalView;
51 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
52 using Element =
typename GridView::template Codim<0>::Entity;
53 using GridIndexType =
typename GridView::IndexSet::IndexType;
54 using Stencil = std::vector<GridIndexType>;
55 using TransmissibilityVector = std::vector<GridIndexType>;
57 using CoordScalar =
typename GridView::ctype;
58 static const int dim = GridView::dimension;
59 static const int dimWorld = GridView::dimensionworld;
61 using FeLocalBasis =
typename GridGeometry::FeCache::FiniteElementType::Traits::LocalBasisType;
62 using ShapeJacobian =
typename FeLocalBasis::Traits::JacobianType;
63 using ShapeValue =
typename Dune::FieldVector<Scalar, 1>;
64 using JacobianInverseTransposed =
typename Element::Geometry::JacobianInverseTransposed;
65 using GlobalPosition = Dune::FieldVector<Scalar, dimWorld>;
70 const Element& element,
71 const FVElementGeometry& fvGeometry,
72 const ElementVolumeVariables& elemVolVars,
73 const SubControlVolumeFace& scvf)
75 const auto geometry = element.geometry();
76 const auto& localBasis = fvGeometry.feLocalBasis();
79 std::vector<ShapeValue> shapeVals;
80 const auto ipLocal = geometry.local(scvf.ipGlobal());
81 jacInvT_ = geometry.jacobianInverseTransposed(ipLocal);
82 localBasis.evaluateJacobian(ipLocal, shapeJacobian_);
83 localBasis.evaluateFunction(ipLocal, shapeVals);
86 shapeValues_.resize(fvGeometry.numScv(), 0.0);
87 if (!scvf.isOnFracture())
88 std::copy(shapeVals.begin(), shapeVals.end(), shapeValues_.begin());
91 const auto thisFacetIdx = scvf.facetIndexInElement();
92 for (
const auto& scv: scvs(fvGeometry))
93 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
94 shapeValues_[scv.indexInElement()] = shapeVals[scv.localDofIndex()];
98 gradN_.resize(fvGeometry.numScv(), GlobalPosition(0.0));
99 if (!scvf.isOnFracture())
101 for (
const auto& scv: scvs(fvGeometry))
102 if (!scv.isOnFracture())
103 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
107 const auto thisFacetIdx = scvf.facetIndexInElement();
110 std::vector<unsigned int> facetLocalDofs;
111 for (
const auto& scv : scvs(fvGeometry))
112 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
113 facetLocalDofs.push_back(scv.localDofIndex());
115 for (
const auto& scv: scvs(fvGeometry))
118 if (scv.isOnFracture() && scv.facetIndexInElement() == thisFacetIdx)
119 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
122 else if (!scv.isOnFracture()
123 && std::find( facetLocalDofs.begin(),
124 facetLocalDofs.end(),
125 scv.localDofIndex() ) == facetLocalDofs.end())
127 jacInvT_.mv(shapeJacobian_[scv.localDofIndex()][0], gradN_[scv.indexInElement()]);
134 const std::vector<ShapeJacobian>&
shapeJacobian()
const {
return shapeJacobian_; }
136 const std::vector<ShapeValue>&
shapeValues()
const {
return shapeValues_; }
138 const JacobianInverseTransposed&
jacInvT()
const {
return jacInvT_; }
140 const GlobalPosition&
gradN(
unsigned int scvIdxInElement)
const {
return gradN_[scvIdxInElement]; }
143 std::vector<GlobalPosition> gradN_;
144 std::vector<ShapeJacobian> shapeJacobian_;
145 std::vector<ShapeValue> shapeValues_;
146 JacobianInverseTransposed jacInvT_;
The available discretization methods in Dumux.
Classes related to flux variables caching.
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
We only store discretization-related quantities for the box method. However, we cannot reuse the cach...
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:43
const GlobalPosition & gradN(unsigned int scvIdxInElement) const
Returns the shape value gradients corresponding to an scv.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:140
const std::vector< ShapeValue > & shapeValues() const
Returns the shape values for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:136
const std::vector< ShapeJacobian > & shapeJacobian() const
Returns the Jacobian of the shape functions at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:134
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:69
const JacobianInverseTransposed & jacInvT() const
Returns the shape value gradients for all scvs at the integration point.
Definition: porousmediumflow/boxdfm/fluxvariablescache.hh:138
Declares all properties used in Dumux.