13#ifndef DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
14#define DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
16#include <dune/common/dynvector.hh>
17#include <dune/common/dynmatrix.hh>
26template<
class TypeTag,
class DiscretizationMethod>
34template<
class TypeTag>
40 using Element =
typename GridView::template Codim<0>::Entity;
42 static constexpr int dim = GridView::dimension;
43 static constexpr int dimWorld = GridView::dimensionworld;
46 using FVElementGeometry =
typename GridGeometry::LocalView;
47 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
51 class MpfaDarcysLawCacheFiller
56 template<
class FluxVariablesCache,
class FluxVariablesCacheFiller>
57 static void fill(FluxVariablesCache& scvfFluxVarsCache,
58 const Problem& problem,
59 const Element& element,
60 const FVElementGeometry& fvGeometry,
61 const ElementVolumeVariables& elemVolVars,
62 const SubControlVolumeFace& scvf,
63 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
66 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
67 scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.secondaryInteractionVolume(),
68 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
69 fluxVarsCacheFiller.secondaryIvDataHandle());
71 scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.primaryInteractionVolume(),
72 fluxVarsCacheFiller.primaryIvLocalFaceData(),
73 fluxVarsCacheFiller.primaryIvDataHandle());
78 class MpfaDarcysLawCache
83 using Stencil =
typename DualGridNodalIndexSet::NodalGridStencilType;
85 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
86 using PrimaryDataHandle =
typename ElementFluxVariablesCache::PrimaryIvDataHandle::AdvectionHandle;
87 using SecondaryDataHandle =
typename ElementFluxVariablesCache::SecondaryIvDataHandle::AdvectionHandle;
90 template<
bool doSecondary = cons
iderSecondaryIVs, std::enable_if_t<doSecondary,
int> = 0 >
91 void setHandlePointer_(
const SecondaryDataHandle& dataHandle)
92 { secondaryHandlePtr_ = &dataHandle; }
95 void setHandlePointer_(
const PrimaryDataHandle& dataHandle)
96 { primaryHandlePtr_ = &dataHandle; }
100 using Filler = MpfaDarcysLawCacheFiller;
110 template<
class IV,
class LocalFaceData,
class DataHandle>
111 void updateAdvection(
const IV& iv,
112 const LocalFaceData& localFaceData,
113 const DataHandle& dataHandle)
115 switchFluxSign_ = localFaceData.isOutsideFace();
116 stencil_ = &iv.stencil();
117 setHandlePointer_(dataHandle.advectionHandle());
121 const Stencil& advectionStencil()
const {
return *stencil_; }
124 const PrimaryDataHandle& advectionPrimaryDataHandle()
const {
return *primaryHandlePtr_; }
125 const SecondaryDataHandle& advectionSecondaryDataHandle()
const {
return *secondaryHandlePtr_; }
128 bool advectionSwitchFluxSign()
const {
return switchFluxSign_; }
131 bool switchFluxSign_;
134 const PrimaryDataHandle* primaryHandlePtr_;
135 const SecondaryDataHandle* secondaryHandlePtr_;
138 const Stencil* stencil_;
159 template<
class ElementFluxVariablesCache>
160 static Scalar
flux(
const Problem& problem,
161 const Element& element,
162 const FVElementGeometry& fvGeometry,
163 const ElementVolumeVariables& elemVolVars,
164 const SubControlVolumeFace& scvf,
165 const unsigned int phaseIdx,
166 const ElementFluxVariablesCache& elemFluxVarsCache)
168 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
171 if (fluxVarsCache.usesSecondaryIv())
172 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionSecondaryDataHandle(), phaseIdx);
174 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionPrimaryDataHandle(), phaseIdx);
178 template<
class Problem,
class FluxVarsCache,
class DataHandle >
179 static Scalar flux_(
const Problem& problem,
180 const FluxVarsCache& cache,
181 const DataHandle& dataHandle,
184 dataHandle.setPhaseIndex(phaseIdx);
186 const bool switchSign = cache.advectionSwitchFluxSign();
188 const auto localFaceIdx = cache.ivLocalFaceIndex();
189 const auto idxInOutside = cache.indexInOutsideFaces();
190 const auto& pj = dataHandle.uj();
191 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
192 : (!switchSign ? dataHandle.T()[localFaceIdx]
193 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
194 Scalar scvfFlux = tij*pj;
197 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
199 scvfFlux += dim == dimWorld ? dataHandle.g()[localFaceIdx]
200 : (!switchSign ? dataHandle.g()[localFaceIdx]
201 : dataHandle.gOutside()[localFaceIdx][idxInOutside]);
static Scalar flux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const unsigned int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the advective flux of a fluid phase across the given sub-control volume face.
Definition: flux/ccmpfa/darcyslaw.hh:160
MpfaDarcysLawCache Cache
Definition: flux/ccmpfa/darcyslaw.hh:147
forward declaration of the method-specific implementation
Definition: flux/ccmpfa/darcyslaw.hh:27
Defines all properties used in Dumux.
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
The available discretization methods in Dumux.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.