Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup CCMpfaFlux | ||
10 | * \brief Darcy's Law for cell-centered finite volume schemes | ||
11 | * with multi-point flux approximation. | ||
12 | */ | ||
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 | |||
19 | #include <dumux/common/properties.hh> | ||
20 | #include <dumux/common/parameters.hh> | ||
21 | #include <dumux/discretization/method.hh> | ||
22 | |||
23 | namespace Dumux { | ||
24 | |||
25 | //! forward declaration of the method-specific implementation | ||
26 | template<class TypeTag, class DiscretizationMethod> | ||
27 | class DarcysLawImplementation; | ||
28 | |||
29 | /*! | ||
30 | * \ingroup CCMpfaFlux | ||
31 | * \brief Darcy's law for cell-centered finite volume schemes | ||
32 | * with multi-point flux approximation. | ||
33 | */ | ||
34 | template<class TypeTag> | ||
35 | class DarcysLawImplementation<TypeTag, DiscretizationMethods::CCMpfa> | ||
36 | { | ||
37 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
38 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
39 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
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 | |||
45 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
46 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
47 | using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace; | ||
48 | using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView; | ||
49 | |||
50 | //! Class that fills the cache corresponding to mpfa Darcy's Law | ||
51 | class MpfaDarcysLawCacheFiller | ||
52 | { | ||
53 | public: | ||
54 | //! Function to fill an MpfaDarcysLawCache of a given scvf | ||
55 | //! This interface has to be met by any advection-related cache filler class | ||
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 |
0/4✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 3 not taken.
|
311491542 | if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex())) |
67 | ✗ | scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.secondaryInteractionVolume(), | |
68 | fluxVarsCacheFiller.secondaryIvLocalFaceData(), | ||
69 | fluxVarsCacheFiller.secondaryIvDataHandle()); | ||
70 | else | ||
71 | 311491542 | scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.primaryInteractionVolume(), | |
72 | fluxVarsCacheFiller.primaryIvLocalFaceData(), | ||
73 | fluxVarsCacheFiller.primaryIvDataHandle()); | ||
74 | ✗ | } | |
75 | }; | ||
76 | |||
77 | //! The cache used in conjunction with the mpfa Darcy's Law | ||
78 | class MpfaDarcysLawCache | ||
79 | { | ||
80 | using DualGridNodalIndexSet = GetPropType<TypeTag, Properties::DualGridNodalIndexSet>; | ||
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 | |||
89 | //! sets the pointer to the data handle (overload for secondary data handles) | ||
90 | template< bool doSecondary = considerSecondaryIVs, std::enable_if_t<doSecondary, int> = 0 > | ||
91 | ✗ | void setHandlePointer_(const SecondaryDataHandle& dataHandle) | |
92 | ✗ | { secondaryHandlePtr_ = &dataHandle; } | |
93 | |||
94 | //! sets the pointer to the data handle (overload for primary data handles) | ||
95 | ✗ | void setHandlePointer_(const PrimaryDataHandle& dataHandle) | |
96 | 311491542 | { primaryHandlePtr_ = &dataHandle; } | |
97 | |||
98 | public: | ||
99 | // export the filler type | ||
100 | using Filler = MpfaDarcysLawCacheFiller; | ||
101 | |||
102 | /*! | ||
103 | * \brief Update cached objects (transmissibilities and gravity). | ||
104 | * This is used for updates with primary interaction volumes. | ||
105 | * | ||
106 | * \param iv The interaction volume this scvf is embedded in | ||
107 | * \param localFaceData iv-local info on this scvf | ||
108 | * \param dataHandle Transmissibility matrix & gravity data of this iv | ||
109 | */ | ||
110 | template<class IV, class LocalFaceData, class DataHandle> | ||
111 | void updateAdvection(const IV& iv, | ||
112 | const LocalFaceData& localFaceData, | ||
113 | const DataHandle& dataHandle) | ||
114 | { | ||
115 | 311491542 | switchFluxSign_ = localFaceData.isOutsideFace(); | |
116 | 311491542 | stencil_ = &iv.stencil(); | |
117 | 622983084 | setHandlePointer_(dataHandle.advectionHandle()); | |
118 | } | ||
119 | |||
120 | //! The stencil corresponding to the transmissibilities (primary type) | ||
121 | ✗ | const Stencil& advectionStencil() const { return *stencil_; } | |
122 | |||
123 | //! The corresponding data handles | ||
124 | ✗ | const PrimaryDataHandle& advectionPrimaryDataHandle() const { return *primaryHandlePtr_; } | |
125 | ✗ | const SecondaryDataHandle& advectionSecondaryDataHandle() const { return *secondaryHandlePtr_; } | |
126 | |||
127 | //! Returns whether or not this scvf is an "outside" face in the scope of the iv. | ||
128 | ✗ | bool advectionSwitchFluxSign() const { return switchFluxSign_; } | |
129 | |||
130 | private: | ||
131 | bool switchFluxSign_; | ||
132 | |||
133 | //! pointers to the corresponding iv-data handles | ||
134 | const PrimaryDataHandle* primaryHandlePtr_; | ||
135 | const SecondaryDataHandle* secondaryHandlePtr_; | ||
136 | |||
137 | //! The stencil, i.e. the grid indices j | ||
138 | const Stencil* stencil_; | ||
139 | }; | ||
140 | |||
141 | public: | ||
142 | using DiscretizationMethod = DiscretizationMethods::CCMpfa; | ||
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 | |||
149 | /*! | ||
150 | * \brief Returns the advective flux of a fluid phase | ||
151 | * across the given sub-control volume face. | ||
152 | * \note This assembles the term | ||
153 | * \f$-|\sigma| \mathbf{n}^T \mathbf{K} \left( \nabla p - \rho \mathbf{g} \right)\f$, | ||
154 | * where \f$|\sigma|\f$ is the area of the face and \f$\mathbf{n}\f$ is the outer | ||
155 | * normal vector. Thus, the flux is given in N*m, and can be converted | ||
156 | * into a volume flux (m^3/s) or mass flux (kg/s) by applying an upwind scheme | ||
157 | * for the mobility or the product of density and mobility, respectively. | ||
158 | */ | ||
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 | |||
177 | private: | ||
178 | template< class Problem, class FluxVarsCache, class DataHandle > | ||
179 | 546950356 | static Scalar flux_(const Problem& problem, | |
180 | const FluxVarsCache& cache, | ||
181 | const DataHandle& dataHandle, | ||
182 | int phaseIdx) | ||
183 | { | ||
184 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 543297476 times.
|
546950356 | dataHandle.setPhaseIndex(phaseIdx); |
185 | |||
186 | 546950356 | const bool switchSign = cache.advectionSwitchFluxSign(); | |
187 | |||
188 | 546950356 | const auto localFaceIdx = cache.ivLocalFaceIndex(); | |
189 | 546950356 | const auto idxInOutside = cache.indexInOutsideFaces(); | |
190 |
2/2✓ Branch 0 taken 141025969 times.
✓ Branch 1 taken 145007949 times.
|
546950356 | const auto& pj = dataHandle.uj(); |
191 |
2/2✓ Branch 0 taken 141025969 times.
✓ Branch 1 taken 145007949 times.
|
807866794 | const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx] |
192 | 423077907 | : (!switchSign ? dataHandle.T()[localFaceIdx] | |
193 | 580031796 | : dataHandle.tijOutside()[localFaceIdx][idxInOutside]); | |
194 | 546950356 | Scalar scvfFlux = tij*pj; | |
195 | |||
196 | // maybe add gravitational acceleration | ||
197 |
6/8✓ Branch 0 taken 143 times.
✓ Branch 1 taken 543297333 times.
✓ Branch 3 taken 41 times.
✓ Branch 4 taken 102 times.
✓ Branch 6 taken 41 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 41 times.
✗ Branch 10 not taken.
|
546950356 | static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity"); |
198 |
2/2✓ Branch 0 taken 399454487 times.
✓ Branch 1 taken 143842989 times.
|
546950356 | if (enableGravity) |
199 |
2/2✓ Branch 0 taken 77109125 times.
✓ Branch 1 taken 79553900 times.
|
1052659076 | scvfFlux += dim == dimWorld ? dataHandle.g()[localFaceIdx] |
200 | 231327375 | : (!switchSign ? dataHandle.g()[localFaceIdx] | |
201 | 318215600 | : dataHandle.gOutside()[localFaceIdx][idxInOutside]); | |
202 | |||
203 | // switch the sign if necessary | ||
204 |
2/2✓ Branch 0 taken 273196513 times.
✓ Branch 1 taken 270100963 times.
|
546950356 | if (switchSign) |
205 | 275012769 | scvfFlux *= -1.0; | |
206 | |||
207 | 546950356 | return scvfFlux; | |
208 | } | ||
209 | }; | ||
210 | |||
211 | } // end namespace | ||
212 | |||
213 | #endif | ||
214 |