24#ifndef DUMUX_DISCRETIZATION_STAGGERED_FREEFLOW_FACEVARIABLES_HH
25#define DUMUX_DISCRETIZATION_STAGGERED_FREEFLOW_FACEVARIABLES_HH
35template<
class Scalar,
int upwindSchemeOrder>
40 std::array<Scalar, upwindSchemeOrder-1>
forward{};
41 std::array<Scalar, upwindSchemeOrder-1>
backward{};
61template<
class FacePrimaryVariables,
int dim,
int upwindSchemeOrder>
64 static constexpr int numPairs = (dim == 2) ? 2 : 4;
65 static constexpr bool useHigherOrder = upwindSchemeOrder > 1;
67 using Scalar =
typename FacePrimaryVariables::block_type;
79 inAxisVelocities_.
self = priVars[0];
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)
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");
103 inAxisVelocities_.
self = faceSol[scvf.dofIndex()];
104 inAxisVelocities_.
opposite = faceSol[scvf.dofIndexOpposingFace()];
106 addHigherOrderInAxisVelocities_(faceSol, scvf, std::integral_constant<bool, useHigherOrder>{});
109 for (
int i = 0; i < velocityParallel_.size(); ++i)
110 velocityParallel_[i].fill(0.0);
112 for (
int i = 0; i < scvf.pairData().size(); ++i)
114 const auto& subFaceData = scvf.pairData(i);
117 velocityLateralInside_[i] = faceSol[subFaceData.lateralPair.first];
119 if (scvf.hasOuterLateral(i))
120 velocityLateralOutside_[i] = faceSol[subFaceData.lateralPair.second];
123 for (
int j = 0; j < upwindSchemeOrder; j++)
125 if (scvf.hasParallelNeighbor(i,j))
126 velocityParallel_[i][j] = faceSol[subFaceData.parallelDofs[j]];
136 return inAxisVelocities_.
self;
152 template<
bool enable = useHigherOrder, std::enable_if_t<enable,
int> = 0>
155 return inAxisVelocities_.
backward[backwardIdx];
163 template<
bool enable = useHigherOrder, std::enable_if_t<enable,
int> = 0>
166 return inAxisVelocities_.
forward[forwardIdx];
177 return velocityParallel_[localSubFaceIdx][parallelDegreeIdx];
187 return velocityLateralInside_[localSubFaceIdx];
197 return velocityLateralOutside_[localSubFaceIdx];
202 template<
class SolVector,
class SubControlVolumeFace>
203 void addHigherOrderInAxisVelocities_(
const SolVector& faceSol,
const SubControlVolumeFace& scvf, std::false_type) {}
205 template<
class SolVector,
class SubControlVolumeFace>
206 void addHigherOrderInAxisVelocities_(
const SolVector& faceSol,
const SubControlVolumeFace& scvf, std::true_type)
211 for (
int i = 0; i < scvf.axisData().inAxisForwardDofs.size(); i++)
213 if (scvf.hasForwardNeighbor(i))
214 inAxisVelocities_.
forward[i]= faceSol[scvf.axisData().inAxisForwardDofs[i]];
219 for (
int i = 0; i < scvf.axisData().inAxisBackwardDofs.size(); i++)
221 if (scvf.hasBackwardNeighbor(i))
222 inAxisVelocities_.
backward[i] = faceSol[scvf.axisData().inAxisBackwardDofs[i]];
226 InAxisVelocities inAxisVelocities_;
227 std::array<std::array<Scalar, upwindSchemeOrder>, numPairs> velocityParallel_;
228 std::array<Scalar, numPairs> velocityLateralInside_;
229 std::array<Scalar, numPairs> velocityLateralOutside_;
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