3.2-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, DiscretizationMethod discMethod>
39class DarcysLawImplementation;
40
46template<class TypeTag>
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:
154 // state the discretization method this implementation belongs to
156
157 // export the type for the corresponding cache
158 using Cache = MpfaDarcysLawCache;
159
161 template<class ElementFluxVariablesCache>
162 static Scalar flux(const Problem& problem,
163 const Element& element,
164 const FVElementGeometry& fvGeometry,
165 const ElementVolumeVariables& elemVolVars,
166 const SubControlVolumeFace& scvf,
167 const unsigned int phaseIdx,
168 const ElementFluxVariablesCache& elemFluxVarsCache)
169 {
170 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
171
172 // forward to the private function taking the iv data handle
173 if (fluxVarsCache.usesSecondaryIv())
174 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionSecondaryDataHandle(), phaseIdx);
175 else
176 return flux_(problem, fluxVarsCache, fluxVarsCache.advectionPrimaryDataHandle(), phaseIdx);
177 }
178
179private:
180 template< class Problem, class FluxVarsCache, class DataHandle >
181 static Scalar flux_(const Problem& problem,
182 const FluxVarsCache& cache,
183 const DataHandle& dataHandle,
184 int phaseIdx)
185 {
186 dataHandle.setPhaseIndex(phaseIdx);
187
188 const bool switchSign = cache.advectionSwitchFluxSign();
189
190 const auto localFaceIdx = cache.ivLocalFaceIndex();
191 const auto idxInOutside = cache.indexInOutsideFaces();
192 const auto& pj = dataHandle.uj();
193 const auto& tij = dim == dimWorld ? dataHandle.T()[localFaceIdx]
194 : (!switchSign ? dataHandle.T()[localFaceIdx]
195 : dataHandle.tijOutside()[localFaceIdx][idxInOutside]);
196 Scalar scvfFlux = tij*pj;
197
198 // maybe add gravitational acceleration
199 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(), "Problem.EnableGravity");
200 if (enableGravity)
201 scvfFlux += dim == dimWorld ? dataHandle.g()[localFaceIdx]
202 : (!switchSign ? dataHandle.g()[localFaceIdx]
203 : dataHandle.gOutside()[localFaceIdx][idxInOutside]);
204
205 // switch the sign if necessary
206 if (switchSign)
207 scvfFlux *= -1.0;
208
209 return scvfFlux;
210 }
211};
212
213} // end namespace
214
215#endif
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
DiscretizationMethod
The available discretization methods in Dumux.
Definition: method.hh:37
Definition: adapt.hh:29
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property (equivalent to old macro GET_PROP_TYPE(....
Definition: propertysystem.hh:149
forward declaration of the method-specific implementation
Definition: flux/darcyslaw.hh:38
MpfaDarcysLawCache Cache
Definition: flux/ccmpfa/darcyslaw.hh:158
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)
Compute the advective flux across an scvf.
Definition: flux/ccmpfa/darcyslaw.hh:162
Declares all properties used in Dumux.