3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
26#define DUMUX_DISCRETIZATION_CC_MPFA_DARCYS_LAW_HH
27
28#include <dune/common/dynvector.hh>
29#include <dune/common/dynmatrix.hh>
30
34
35namespace Dumux {
36
38template<class TypeTag, class DiscretizationMethod>
39class DarcysLawImplementation;
40
46template<class TypeTag>
47class DarcysLawImplementation<TypeTag, DiscretizationMethods::CCMpfa>
48{
52 using Element = typename GridView::template Codim<0>::Entity;
53
54 static constexpr int dim = GridView::dimension;
55 static constexpr int dimWorld = GridView::dimensionworld;
56
58 using FVElementGeometry = typename GridGeometry::LocalView;
59 using SubControlVolumeFace = typename FVElementGeometry::SubControlVolumeFace;
60 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
61
63 class MpfaDarcysLawCacheFiller
64 {
65 public:
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)
76 {
77 // get interaction volume related data from the filler class & upate the cache
78 if (fvGeometry.gridGeometry().vertexUsesSecondaryInteractionVolume(scvf.vertexIndex()))
79 scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.secondaryInteractionVolume(),
80 fluxVarsCacheFiller.secondaryIvLocalFaceData(),
81 fluxVarsCacheFiller.secondaryIvDataHandle());
82 else
83 scvfFluxVarsCache.updateAdvection(fluxVarsCacheFiller.primaryInteractionVolume(),
84 fluxVarsCacheFiller.primaryIvLocalFaceData(),
85 fluxVarsCacheFiller.primaryIvDataHandle());
86 }
87 };
88
90 class MpfaDarcysLawCache
91 {
93 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
94
95 using Stencil = typename DualGridNodalIndexSet::NodalGridStencilType;
96
97 static constexpr bool considerSecondaryIVs = GridGeometry::MpfaHelper::considerSecondaryIVs();
98 using PrimaryDataHandle = typename ElementFluxVariablesCache::PrimaryIvDataHandle::AdvectionHandle;
99 using SecondaryDataHandle = typename ElementFluxVariablesCache::SecondaryIvDataHandle::AdvectionHandle;
100
102 template< bool doSecondary = considerSecondaryIVs, std::enable_if_t<doSecondary, int> = 0 >
103 void setHandlePointer_(const SecondaryDataHandle& dataHandle)
104 { secondaryHandlePtr_ = &dataHandle; }
105
107 void setHandlePointer_(const PrimaryDataHandle& dataHandle)
108 { primaryHandlePtr_ = &dataHandle; }
109
110 public:
111 // export the filler type
112 using Filler = MpfaDarcysLawCacheFiller;
113
122 template<class IV, class LocalFaceData, class DataHandle>
123 void updateAdvection(const IV& iv,
124 const LocalFaceData& localFaceData,
125 const DataHandle& dataHandle)
126 {
127 switchFluxSign_ = localFaceData.isOutsideFace();
128 stencil_ = &iv.stencil();
129 setHandlePointer_(dataHandle.advectionHandle());
130 }
131
133 const Stencil& advectionStencil() const { return *stencil_; }
134
136 const PrimaryDataHandle& advectionPrimaryDataHandle() const { return *primaryHandlePtr_; }
137 const SecondaryDataHandle& advectionSecondaryDataHandle() const { return *secondaryHandlePtr_; }
138
140 bool advectionSwitchFluxSign() const { return switchFluxSign_; }
141
142 private:
143 bool switchFluxSign_;
144
146 const PrimaryDataHandle* primaryHandlePtr_;
147 const SecondaryDataHandle* secondaryHandlePtr_;
148
150 const Stencil* stencil_;
151 };
152
153public:
155 // state the discretization method this implementation belongs to
156 static constexpr DiscretizationMethod discMethod{};
157
158 // export the type for the corresponding cache
159 using Cache = MpfaDarcysLawCache;
160
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)
179 {
180 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
181
182 // forward to the private function taking the iv data handle
183 if (fluxVarsCache.usesSecondaryIv())
184 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionSecondaryDataHandle(), phaseIdx);
185 else
186 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionPrimaryDataHandle(), phaseIdx);
187 }
188
189private:
190 template< class Problem, class FluxVarsCache, class DataHandle >
191 static Scalar flux_(const Problem& problem,
192 const FluxVarsCache& cache,
193 const DataHandle& dataHandle,
194 int phaseIdx)
195 {
196 dataHandle.setPhaseIndex(phaseIdx);
197
198 const bool switchSign = cache.advectionSwitchFluxSign();
199
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;
207
208 // maybe add gravitational acceleration
209 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
210 if (enableGravity)
211 scvfFlux += dim == dimWorld ? dataHandle.g()[localFaceIdx]
212 : (!switchSign ? dataHandle.g()[localFaceIdx]
213 : dataHandle.gOutside()[localFaceIdx][idxInOutside]);
214
215 // switch the sign if necessary
216 if (switchSign)
217 scvfFlux *= -1.0;
218
219 return scvfFlux;
220 }
221};
222
223} // end namespace
224
225#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
Definition: method.hh:61
forward declaration of the method-specific implementation
Definition: flux/box/darcyslaw.hh:42
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.