25#ifndef DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
26#define DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
28#include <dune/common/dynvector.hh>
29#include <dune/common/dynmatrix.hh>
38template<
class TypeTag,
class DiscretizationMethod>
39class DarcysLawImplementation;
46template<
class TypeTag>
52 using Element =
typename GridView::template Codim<0>::Entity;
54 static constexpr int dim = GridView::dimension;
55 static constexpr int dimWorld = GridView::dimensionworld;
58 using FVElementGeometry =
typename GridGeometry::LocalView;
59 using SubControlVolumeFace =
typename FVElementGeometry::SubControlVolumeFace;
63 class MpfaDarcysLawCacheFiller
68 template<
class FluxVariablesCache,
class FluxVariablesCacheFiller>
69 static void fill(FluxVariablesCache& scvfFluxVarsCache,
70 const Problem& problem,
71 const Element& element,
72 const FVElementGeometry& fvGeometry,
73 const ElementVolumeVariables& elemVolVars,
74 const SubControlVolumeFace& scvf,
75 const FluxVariablesCacheFiller& fluxVarsCacheFiller)
78 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
79 scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.secondaryInteractionVolume(),
80 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
81 fluxVarsCacheFiller.secondaryIvDataHandle());
83 scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.primaryInteractionVolume(),
84 fluxVarsCacheFiller.primaryIvLocalFaceData(),
85 fluxVarsCacheFiller.primaryIvDataHandle());
90 class MpfaDarcysLawCache
95 using Stencil =
typename DualGridNodalIndexSet::NodalGridStencilType;
97 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
98 using PrimaryDataHandle =
typename ElementFluxVariablesCache::PrimaryIvDataHandle::AdvectionHandle;
99 using SecondaryDataHandle =
typename ElementFluxVariablesCache::SecondaryIvDataHandle::AdvectionHandle;
102 template<
bool doSecondary = cons
iderSecondaryIVs, std::enable_if_t<doSecondary,
int> = 0 >
103 void setHandlePointer_(
const SecondaryDataHandle& dataHandle)
104 { secondaryHandlePtr_ = &dataHandle; }
107 void setHandlePointer_(
const PrimaryDataHandle& dataHandle)
108 { primaryHandlePtr_ = &dataHandle; }
112 using Filler = MpfaDarcysLawCacheFiller;
122 template<
class IV,
class LocalFaceData,
class DataHandle>
123 void updateAdvection(
const IV& iv,
124 const LocalFaceData& localFaceData,
125 const DataHandle& dataHandle)
127 switchFluxSign_ = localFaceData.isOutsideFace();
128 stencil_ = &iv.stencil();
129 setHandlePointer_(dataHandle.advectionHandle());
133 const Stencil& advectionStencil()
const {
return *stencil_; }
136 const PrimaryDataHandle& advectionPrimaryDataHandle()
const {
return *primaryHandlePtr_; }
137 const SecondaryDataHandle& advectionSecondaryDataHandle()
const {
return *secondaryHandlePtr_; }
140 bool advectionSwitchFluxSign()
const {
return switchFluxSign_; }
143 bool switchFluxSign_;
146 const PrimaryDataHandle* primaryHandlePtr_;
147 const SecondaryDataHandle* secondaryHandlePtr_;
150 const Stencil* stencil_;
171 template<
class ElementFluxVariablesCache>
172 static Scalar
flux(
const Problem& problem,
173 const Element& element,
174 const FVElementGeometry& fvGeometry,
175 const ElementVolumeVariables& elemVolVars,
176 const SubControlVolumeFace& scvf,
177 const unsigned int phaseIdx,
178 const ElementFluxVariablesCache& elemFluxVarsCache)
180 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
183 if (fluxVarsCache.usesSecondaryIv())
184 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionSecondaryDataHandle(), phaseIdx);
186 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionPrimaryDataHandle(), phaseIdx);
190 template<
class Problem,
class FluxVarsCache,
class DataHandle >
191 static Scalar flux_(
const Problem& problem,
192 const FluxVarsCache& cache,
193 const DataHandle& dataHandle,
196 dataHandle.setPhaseIndex(phaseIdx);
198 const bool switchSign = cache.advectionSwitchFluxSign();
200 const auto localFaceIdx = cache.ivLocalFaceIndex();
201 const auto idxInOutside = cache.indexInOutsideFaces();
202 const auto& pj = dataHandle.uj();
203 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
204 : (!switchSign ? dataHandle.T()[localFaceIdx]
205 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
206 Scalar scvfFlux = tij*pj;
209 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
211 scvfFlux += dim == dimWorld ? dataHandle.g()[localFaceIdx]
212 : (!switchSign ? dataHandle.g()[localFaceIdx]
213 : dataHandle.gOutside()[localFaceIdx][idxInOutside]);
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:44
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:172
MpfaDarcysLawCache Cache
Definition: flux/ccmpfa/darcyslaw.hh:159
Declares all properties used in Dumux.