3.2-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#include <type_traits>
30
31namespace Dumux {
32
33namespace Detail {
34
35template<class Scalar, int upwindSchemeOrder>
37{
38 Scalar self = 0.0;
39 Scalar opposite = 0.0;
40 std::array<Scalar, upwindSchemeOrder-1> forward{};
41 std::array<Scalar, upwindSchemeOrder-1> backward{};
42};
43
44template<class Scalar>
45struct InAxisVelocities<Scalar, 1>
46{
47 Scalar self = 0.0;
48 Scalar opposite = 0.0;
49};
50
51} // end namespace Detail
52
61template<class FacePrimaryVariables, int dim, int upwindSchemeOrder>
63{
64 static constexpr int numPairs = (dim == 2) ? 2 : 4;
65 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
66
67 using Scalar = typename FacePrimaryVariables::block_type;
69
70public:
71
77 void updateOwnFaceOnly(const FacePrimaryVariables& priVars)
78 {
79 inAxisVelocities_.self = priVars[0];
80 }
81
92 template<class FaceSolution, class Problem, class Element,
93 class FVElementGeometry, class SubControlVolumeFace>
94 void update(const FaceSolution& faceSol,
95 const Problem& problem,
96 const Element& element,
97 const FVElementGeometry& fvGeometry,
98 const SubControlVolumeFace& scvf)
99 {
100 static_assert(std::decay_t<decltype(faceSol[0])>::dimension == 1,
101 "\n\n\nVelocity primary variable must be a scalar value. \n\n Make sure to use\n\n ffSol = partial(sol, ffFaceIdx, ffCellCenterIdx);\n\n");
102
103 inAxisVelocities_.self = faceSol[scvf.dofIndex()];
104 inAxisVelocities_.opposite = faceSol[scvf.dofIndexOpposingFace()];
105
106 addHigherOrderInAxisVelocities_(faceSol, scvf, std::integral_constant<bool, useHigherOrder>{});
107
108 // handle all sub faces
109 for (int i = 0; i < velocityParallel_.size(); ++i)
110 velocityParallel_[i].fill(0.0);
111
112 for (int i = 0; i < scvf.pairData().size(); ++i)
113 {
114 const auto& subFaceData = scvf.pairData(i);
115
116 // treat the velocities normal to the face
117 velocityLateralInside_[i] = faceSol[subFaceData.lateralPair.first];
118
119 if (scvf.hasOuterLateral(i))
120 velocityLateralOutside_[i] = faceSol[subFaceData.lateralPair.second];
121
122 // treat the velocities parallel to the self face
123 for (int j = 0; j < upwindSchemeOrder; j++)
124 {
125 if (scvf.hasParallelNeighbor(i,j))
126 velocityParallel_[i][j] = faceSol[subFaceData.parallelDofs[j]];
127 }
128 }
129 }
130
134 Scalar velocitySelf() const
135 {
136 return inAxisVelocities_.self;
137 }
138
142 Scalar velocityOpposite() const
143 {
144 return inAxisVelocities_.opposite;
145 }
146
152 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
153 Scalar velocityBackward(const int backwardIdx) const
154 {
155 return inAxisVelocities_.backward[backwardIdx];
156 }
157
163 template<bool enable = useHigherOrder, std::enable_if_t<enable, int> = 0>
164 Scalar velocityForward(const int forwardIdx) const
165 {
166 return inAxisVelocities_.forward[forwardIdx];
167 }
168
175 Scalar velocityParallel(const int localSubFaceIdx, const int parallelDegreeIdx = 0) const
176 {
177 return velocityParallel_[localSubFaceIdx][parallelDegreeIdx];
178 }
179
185 Scalar velocityLateralInside(const int localSubFaceIdx) const
186 {
187 return velocityLateralInside_[localSubFaceIdx];
188 }
189
195 Scalar velocityLateralOutside(const int localSubFaceIdx) const
196 {
197 return velocityLateralOutside_[localSubFaceIdx];
198 }
199
200private:
201
202 template<class SolVector, class SubControlVolumeFace>
203 void addHigherOrderInAxisVelocities_(const SolVector& faceSol, const SubControlVolumeFace& scvf, std::false_type) {}
204
205 template<class SolVector, class SubControlVolumeFace>
206 void addHigherOrderInAxisVelocities_(const SolVector& faceSol, const SubControlVolumeFace& scvf, std::true_type)
207 {
208
209 // treat the velocity forward of the self face i.e. the face that is
210 // forward wrt the self face by degree i
211 for (int i = 0; i < scvf.axisData().inAxisForwardDofs.size(); i++)
212 {
213 if (scvf.hasForwardNeighbor(i))
214 inAxisVelocities_.forward[i]= faceSol[scvf.axisData().inAxisForwardDofs[i]];
215 }
216
217 // treat the velocity at the first backward face i.e. the face that is
218 // behind the opposite face by degree i
219 for (int i = 0; i < scvf.axisData().inAxisBackwardDofs.size(); i++)
220 {
221 if (scvf.hasBackwardNeighbor(i))
222 inAxisVelocities_.backward[i] = faceSol[scvf.axisData().inAxisBackwardDofs[i]];
223 }
224 }
225
226 InAxisVelocities inAxisVelocities_;
227 std::array<std::array<Scalar, upwindSchemeOrder>, numPairs> velocityParallel_;
228 std::array<Scalar, numPairs> velocityLateralInside_;
229 std::array<Scalar, numPairs> velocityLateralOutside_;
230
231};
232
233} // end namespace Dumux
234
235#endif
Definition: adapt.hh:29
Definition: facevariables.hh:37
std::array< Scalar, upwindSchemeOrder-1 > forward
Definition: facevariables.hh:40
Scalar opposite
Definition: facevariables.hh:39
std::array< Scalar, upwindSchemeOrder-1 > backward
Definition: facevariables.hh:41
Scalar self
Definition: facevariables.hh:38
The face variables class for free flow staggered grid models. Contains all relevant velocities for th...
Definition: facevariables.hh:63
Scalar velocityForward(const int forwardIdx) const
Returns the velocity at a forward face.
Definition: facevariables.hh:164
void update(const FaceSolution &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:94
Scalar velocityOpposite() const
Returns the velocity at the opposing face.
Definition: facevariables.hh:142
Scalar velocityParallel(const int localSubFaceIdx, const int parallelDegreeIdx=0) const
Returns the velocity at a parallel face.
Definition: facevariables.hh:175
void updateOwnFaceOnly(const FacePrimaryVariables &priVars)
Partial update of the face variables. Only the face itself is considered.
Definition: facevariables.hh:77
Scalar velocitySelf() const
Returns the velocity at the face itself.
Definition: facevariables.hh:134
Scalar velocityLateralInside(const int localSubFaceIdx) const
Returns the velocity at the inner normal face.
Definition: facevariables.hh:185
Scalar velocityLateralOutside(const int localSubFaceIdx) const
Returns the velocity at the outer normal face.
Definition: facevariables.hh:195
Scalar velocityBackward(const int backwardIdx) const
Returns the velocity at a backward face.
Definition: facevariables.hh:153