version 3.8
flux/ccmpfa/darcyslaw.hh
Go to the documentation of this file.
1// -*- mode: C++; tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 4 -*-
2// vi: set et ts=4 sw=4 sts=4:
3//
4// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
14#define DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
15
16#include <dune/common/dynvector.hh>
17#include <dune/common/dynmatrix.hh>
18
22
23namespace Dumux {
24
26template<class TypeTag, class DiscretizationMethod>
28
34template<class TypeTag>
35class DarcysLawImplementation<TypeTag, DiscretizationMethods::CCMpfa>
36{
40 using Element = typename GridView::template Codim<0>::Entity;
41
42 static constexpr int dim = GridView::dimension;
43 static constexpr int dimWorld = GridView::dimensionworld;
44
46 using FVElementGeometry = typename GridGeometry::LocalView;
47 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
48 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
49
51 class MpfaDarcysLawCacheFiller
52 {
53 public:
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)
64 {
65 // get interaction volume related data from the filler class & update the cache
66 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
67 scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.secondaryInteractionVolume(),
68 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
69 fluxVarsCacheFiller.secondaryIvDataHandle());
70 else
71 scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.primaryInteractionVolume(),
72 fluxVarsCacheFiller.primaryIvLocalFaceData(),
73 fluxVarsCacheFiller.primaryIvDataHandle());
74 }
75 };
76
78 class MpfaDarcysLawCache
79 {
81 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
82
83 using Stencil = typename DualGridNodalIndexSet::NodalGridStencilType;
84
85 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
86 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle::AdvectionHandle;
87 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle::AdvectionHandle;
88
90 template< bool doSecondary = considerSecondaryIVs, std::enable_if_t<doSecondary, int> = 0 >
91 void setHandlePointer_(const SecondaryDataHandle& dataHandle)
92 { secondaryHandlePtr_ = &dataHandle; }
93
95 void setHandlePointer_(const PrimaryDataHandle& dataHandle)
96 { primaryHandlePtr_ = &dataHandle; }
97
98 public:
99 // export the filler type
100 using Filler = MpfaDarcysLawCacheFiller;
101
110 template<class IV, class LocalFaceData, class DataHandle>
111 void updateAdvection(const IV& iv,
112 const LocalFaceData& localFaceData,
113 const DataHandle& dataHandle)
114 {
115 switchFluxSign_ = localFaceData.isOutsideFace();
116 stencil_ = &iv.stencil();
117 setHandlePointer_(dataHandle.advectionHandle());
118 }
119
121 const Stencil& advectionStencil() const { return *stencil_; }
122
124 const PrimaryDataHandle& advectionPrimaryDataHandle() const { return *primaryHandlePtr_; }
125 const SecondaryDataHandle& advectionSecondaryDataHandle() const { return *secondaryHandlePtr_; }
126
128 bool advectionSwitchFluxSign() const { return switchFluxSign_; }
129
130 private:
131 bool switchFluxSign_;
132
134 const PrimaryDataHandle* primaryHandlePtr_;
135 const SecondaryDataHandle* secondaryHandlePtr_;
136
138 const Stencil* stencil_;
139 };
140
141public:
143 // state the discretization method this implementation belongs to
144 static constexpr DiscretizationMethod discMethod{};
145
146 // export the type for the corresponding cache
147 using Cache = MpfaDarcysLawCache;
148
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)
167 {
168 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
169
170 // forward to the private function taking the iv data handle
171 if (fluxVarsCache.usesSecondaryIv())
172 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionSecondaryDataHandle(), phaseIdx);
173 else
174 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionPrimaryDataHandle(), phaseIdx);
175 }
176
177private:
178 template< class Problem, class FluxVarsCache, class DataHandle >
179 static Scalar flux_(const Problem& problem,
180 const FluxVarsCache& cache,
181 const DataHandle& dataHandle,
182 int phaseIdx)
183 {
184 dataHandle.setPhaseIndex(phaseIdx);
185
186 const bool switchSign = cache.advectionSwitchFluxSign();
187
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;
195
196 // maybe add gravitational acceleration
197 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
198 if (enableGravity)
199 scvfFlux += dim == dimWorld ? dataHandle.g()[localFaceIdx]
200 : (!switchSign ? dataHandle.g()[localFaceIdx]
201 : dataHandle.gOutside()[localFaceIdx][idxInOutside]);
202
203 // switch the sign if necessary
204 if (switchSign)
205 scvfFlux *= -1.0;
206
207 return scvfFlux;
208 }
209};
210
211} // end namespace
212
213#endif
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.
Definition: adapt.hh:17
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Definition: method.hh:33