3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
facevariables.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 *****************************************************************************/
24#ifndef DUMUX_DISCRETIZATION_STAGGERED_FREEFLOW_FACEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FREEFLOW_FACEVARIABLES_HH
26
27#include <array>
28#include <vector>
29
30namespace Dumux {
31
32namespace Detail {
33
34template<class Scalar, int upwindSchemeOrder>
36{
37 Scalar self = 0.0;
38 Scalar opposite = 0.0;
39 std::array<Scalar, upwindSchemeOrder-1> forward{};
40 std::array<Scalar, upwindSchemeOrder-1> backward{};
41};
42
43template<class Scalar>
44struct InAxisVelocities<Scalar, 1>
45{
46 Scalar self = 0.0;
47 Scalar opposite = 0.0;
48};
49
50} // end namespace Detail
51
60template<class FacePrimaryVariables, int dim, int upwindSchemeOrder>
62{
63 static constexpr int numPairs = (dim == 2) ? 2 : 4;
64 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
65
66 using Scalar = typename FacePrimaryVariables::block_type;
68
69public:
70
76 void updateOwnFaceOnly(const FacePrimaryVariables& priVars)
77 {
78 inAxisVelocities_.self = priVars[0];
79 }
80
91 template<class SolVector, class Problem, class Element,
92 class FVElementGeometry, class SubControlVolumeFace>
93 void update(const SolVector& faceSol,
94 const Problem& problem,
95 const Element& element,
96 const FVElementGeometry& fvGeometry,
97 const SubControlVolumeFace& scvf)
98 {
99 inAxisVelocities_.self = faceSol[scvf.dofIndex()];
100 inAxisVelocities_.opposite = faceSol[scvf.dofIndexOpposingFace()];
101
102 addHigherOrderInAxisVelocities_(faceSol, scvf, std::integral_constant<bool, useHigherOrder>{});
103
104 // handle all sub faces
105 for (int i = 0; i < velocityParallel_.size(); ++i)
106 velocityParallel_[i].fill(0.0);
107
108 for (int i = 0; i < scvf.pairData().size(); ++i)
109 {
110 const auto& subFaceData = scvf.pairData(i);
111
112 // treat the velocities normal to the face
113 velocityLateralInside_[i] = faceSol[subFaceData.lateralPair.first];
114
115 if (scvf.hasOuterLateral(i))
116 velocityLateralOutside_[i] = faceSol[subFaceData.lateralPair.second];
117
118 // treat the velocities parallel to the self face
119 for (int j = 0; j < upwindSchemeOrder; j++)
120 {
121 if (scvf.hasParallelNeighbor(i,j))
122 velocityParallel_[i][j] = faceSol[subFaceData.parallelDofs[j]];
123 }
124 }
125 }
126
130 Scalar velocitySelf() const
131 {
132 return inAxisVelocities_.self;
133 }
134
138 Scalar velocityOpposite() const
139 {
140 return inAxisVelocities_.opposite;
141 }
142
148 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
149 Scalar velocityBackward(const int backwardIdx) const
150 {
151 return inAxisVelocities_.backward[backwardIdx];
152 }
153
159 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
160 Scalar velocityForward(const int forwardIdx) const
161 {
162 return inAxisVelocities_.forward[forwardIdx];
163 }
164
171 Scalar velocityParallel(const int localSubFaceIdx, const int parallelDegreeIdx = 0) const
172 {
173 return velocityParallel_[localSubFaceIdx][parallelDegreeIdx];
174 }
175
181 Scalar velocityLateralInside(const int localSubFaceIdx) const
182 {
183 return velocityLateralInside_[localSubFaceIdx];
184 }
185
191 Scalar velocityLateralOutside(const int localSubFaceIdx) const
192 {
193 return velocityLateralOutside_[localSubFaceIdx];
194 }
195
196private:
197
198 template<class SolVector, class SubControlVolumeFace>
199 void addHigherOrderInAxisVelocities_(const SolVector& faceSol, const SubControlVolumeFace& scvf, std::false_type) {}
200
201 template<class SolVector, class SubControlVolumeFace>
202 void addHigherOrderInAxisVelocities_(const SolVector& faceSol, const SubControlVolumeFace& scvf, std::true_type)
203 {
204
205 // treat the velocity forward of the self face i.e. the face that is
206 // forward wrt the self face by degree i
207 for (int i = 0; i < scvf.axisData().inAxisForwardDofs.size(); i++)
208 {
209 if (scvf.hasForwardNeighbor(i))
210 inAxisVelocities_.forward[i]= faceSol[scvf.axisData().inAxisForwardDofs[i]];
211 }
212
213 // treat the velocity at the first backward face i.e. the face that is
214 // behind the opposite face by degree i
215 for (int i = 0; i < scvf.axisData().inAxisBackwardDofs.size(); i++)
216 {
217 if (scvf.hasBackwardNeighbor(i))
218 inAxisVelocities_.backward[i] = faceSol[scvf.axisData().inAxisBackwardDofs[i]];
219 }
220 }
221
222 InAxisVelocities inAxisVelocities_;
223 std::array<std::array<Scalar, upwindSchemeOrder>, numPairs> velocityParallel_;
224 std::array<Scalar, numPairs> velocityLateralInside_;
225 std::array<Scalar, numPairs> velocityLateralOutside_;
226
227};
228
229} // end namespace Dumux
230
231#endif
make the local view function available whenever we use the grid geometry
Definition: adapt.hh:29
Definition: facevariables.hh:36
std::array< Scalar, upwindSchemeOrder-1 > forward
Definition: facevariables.hh:39
Scalar opposite
Definition: facevariables.hh:38
std::array< Scalar, upwindSchemeOrder-1 > backward
Definition: facevariables.hh:40
Scalar self
Definition: facevariables.hh:37
The face variables class for free flow staggered grid models. Contains all relevant velocities for th...
Definition: facevariables.hh:62
Scalar velocityForward(const int forwardIdx) const
Returns the velocity at a forward face.
Definition: facevariables.hh:160
void update(const SolVector &faceSol, const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const SubControlVolumeFace &scvf)
Complete update of the face variables (i.e. velocities for free flow) for a given face.
Definition: facevariables.hh:93
Scalar velocityOpposite() const
Returns the velocity at the opposing face.
Definition: facevariables.hh:138
Scalar velocityParallel(const int localSubFaceIdx, const int parallelDegreeIdx=0) const
Returns the velocity at a parallel face.
Definition: facevariables.hh:171
void updateOwnFaceOnly(const FacePrimaryVariables &priVars)
Partial update of the face variables. Only the face itself is considered.
Definition: facevariables.hh:76
Scalar velocitySelf() const
Returns the velocity at the face itself.
Definition: facevariables.hh:130
Scalar velocityLateralInside(const int localSubFaceIdx) const
Returns the velocity at the inner normal face.
Definition: facevariables.hh:181
Scalar velocityLateralOutside(const int localSubFaceIdx) const
Returns the velocity at the outer normal face.
Definition: facevariables.hh:191
Scalar velocityBackward(const int backwardIdx) const
Returns the velocity at a backward face.
Definition: facevariables.hh:149