26#ifndef DUMUX_FLUX_FORCHHEIMER_VELOCITY_HH
27#define DUMUX_FLUX_FORCHHEIMER_VELOCITY_HH
29#include <dune/common/fvector.hh>
30#include <dune/common/fmatrix.hh>
47template<
class Scalar,
class Gr
idGeometry,
class FluxVariables>
50 using FVElementGeometry =
typename GridGeometry::LocalView;
51 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
52 using GridView =
typename GridGeometry::GridView;
53 using Element =
typename GridView::template Codim<0>::Entity;
55 static constexpr int dim = GridView::dimension;
56 static constexpr int dimWorld = GridView::dimensionworld;
57 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
122 template<
class ElementVolumeVariables>
124 const ElementVolumeVariables& elemVolVars,
125 const SubControlVolumeFace& scvf,
130 const auto& problem = elemVolVars.gridVolVars().problem();
133 const Scalar forchCoeff = problem.spatialParams().forchCoeff(scvf);
143 static const Scalar epsilon = getParamFromGroup<Scalar>(problem.paramGroup(),
"Forchheimer.NewtonTolerance", 1e-12);
144 static const std::size_t maxNumIter = getParamFromGroup<std::size_t>(problem.paramGroup(),
"Forchheimer.MaxIterations", 30);
145 for (
int k = 0; residual.two_norm() > epsilon ; ++k)
148 DUNE_THROW(
NumericalProblem,
"could not determine forchheimer velocity within "
149 << k <<
" iterations");
152 forchheimerResidual_(residual,
162 forchheimerDerivative_(gradF,
171 gradF.solve(deltaV, residual);
177 auto upwindTerm = [phaseIdx](
const auto& volVars){
return volVars.mobility(phaseIdx); };
178 const Scalar forchheimerUpwindMobility = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm,
velocity * scvf.unitOuterNormal(), phaseIdx);
181 if (forchheimerUpwindMobility > 0.0)
182 velocity /= forchheimerUpwindMobility;
191 template<
class Problem,
class ElementVolumeVariables,
192 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
193 std::enable_if_t<scalarPerm, int> = 0>
195 const ElementVolumeVariables& elemVolVars,
196 const SubControlVolumeFace& scvf)
199 Scalar harmonicMeanSqrtK(0.0);
201 const auto insideScvIdx = scvf.insideScvIdx();
202 const auto& insideVolVars = elemVolVars[insideScvIdx];
203 const Scalar Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
204 const Scalar sqrtKi = sqrt(Ki);
206 if (!scvf.boundary())
208 const auto outsideScvIdx = scvf.outsideScvIdx();
209 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
210 const Scalar Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
211 const Scalar sqrtKj = sqrt(Kj);
215 harmonicMeanSqrtK = sqrtKi;
219 for (
int i = 0; i < dimWorld; ++i)
220 result[i][i] = harmonicMeanSqrtK;
229 template<
class Problem,
class ElementVolumeVariables,
230 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
231 std::enable_if_t<!scalarPerm, int> = 0>
233 const ElementVolumeVariables& elemVolVars,
234 const SubControlVolumeFace& scvf)
241 const auto insideScvIdx = scvf.insideScvIdx();
242 const auto& insideVolVars = elemVolVars[insideScvIdx];
243 const auto& Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
246 if (!isDiagonal_(Ki))
247 DUNE_THROW(Dune::NotImplemented,
"Only diagonal permeability tensors are supported.");
249 for (
int i = 0; i < dim; ++i)
250 sqrtKi[i][i] = sqrt(Ki[i][i]);
252 if (!scvf.boundary())
254 const auto outsideScvIdx = scvf.outsideScvIdx();
255 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
256 const auto& Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
258 if (!isDiagonal_(Kj))
259 DUNE_THROW(Dune::NotImplemented,
"Only diagonal permeability tensors are supported.");
261 for (
int i = 0; i < dim; ++i)
262 sqrtKj[i][i] = sqrt(Kj[i][i]);
267 harmonicMeanSqrtK = sqrtKi;
269 return harmonicMeanSqrtK;
275 template<
class ElementVolumeVariables>
277 const Scalar forchCoeff,
281 const ElementVolumeVariables& elemVolVars,
282 const SubControlVolumeFace& scvf,
285 residual = oldVelocity;
287 residual -= darcyVelocity;
289 auto upwindTerm = [phaseIdx](
const auto& volVars){
return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
290 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, oldVelocity * scvf.unitOuterNormal(), phaseIdx);
291 sqrtK.usmv(forchCoeff * rhoOverMu * oldVelocity.two_norm(), oldVelocity, residual);
295 template<
class ElementVolumeVariables>
297 const Scalar forchCoeff,
300 const ElementVolumeVariables& elemVolVars,
301 const SubControlVolumeFace& scvf,
307 auto upwindTerm = [phaseIdx](
const auto& volVars){
return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
309 const Scalar absV =
velocity.two_norm();
319 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm,
velocity * scvf.unitOuterNormal(), phaseIdx);
322 for (
int i = 0; i < dim; i++)
324 for (
int k = 0; k <= i; k++)
330 derivative[k][i] = derivative[i][k];
337 for (
int i = 0; i < dim; i++)
338 derivative[i][i] += (1.0 + (sqrtK[i][i]*forchCoeff * absV * rhoOverMu));
352 for (
int i = 0; i < dim; i++)
353 for (
int k = 0; k < dim; k++)
354 if ((i != k) && (abs(K[i][k]) >= 1e-25))
359 template<
class Problem,
class VolumeVariables>
360 static decltype(
auto) getPermeability_(
const Problem& problem,
361 const VolumeVariables& volVars,
362 const GlobalPosition& scvfIpGlobal)
364 if constexpr (Problem::SpatialParams::evaluatePermeabilityAtScvfIP())
365 return problem.spatialParams().permeabilityAtPos(scvfIpGlobal);
367 return volVars.permeability();
Define some often used mathematical functions.
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
A free function to average a Tensor at an interface.
UpwindSchemeImpl< GridGeometry, typename GridGeometry::DiscretizationMethod > UpwindScheme
The upwind scheme used for the advective fluxes. This depends on the chosen discretization method.
Definition: flux/upwindscheme.hh:42
Scalar faceTensorAverage(const Scalar T1, const Scalar T2, const Dune::FieldVector< Scalar, dim > &normal)
Average of a discontinuous scalar field at discontinuity interface (for compatibility reasons with th...
Definition: facetensoraverage.hh:41
Exception thrown if a fixable numerical problem occurs.
Definition: exceptions.hh:39
Forchheimer's law This file contains the calculation of the Forchheimer velocity for a given Darcy ve...
Definition: forchheimervelocity.hh:49
typename FluxVariables::UpwindScheme UpwindScheme
Definition: forchheimervelocity.hh:59
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: forchheimervelocity.hh:194
Dune::FieldVector< Scalar, dimWorld > DimWorldVector
Definition: forchheimervelocity.hh:60
static DimWorldVector velocity(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const DimWorldMatrix sqrtK, const DimWorldVector darcyVelocity)
Compute the Forchheimer velocity of a phase at the given sub-control volume face using the Forchheime...
Definition: forchheimervelocity.hh:123
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimWorldMatrix
Definition: forchheimervelocity.hh:61
Declares all properties used in Dumux.