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>
49template<
class Scalar,
class Gr
idGeometry,
class FluxVariables>
52 using FVElementGeometry =
typename GridGeometry::LocalView;
53 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
54 using GridView =
typename GridGeometry::GridView;
55 using Element =
typename GridView::template Codim<0>::Entity;
57 static constexpr int dim = GridView::dimension;
58 static constexpr int dimWorld = GridView::dimensionworld;
59 using GlobalPosition =
typename Element::Geometry::GlobalCoordinate;
66 template<
class ElementVolumeVariables>
68 const ElementVolumeVariables& elemVolVars,
69 const SubControlVolumeFace& scvf,
74 const auto& problem = elemVolVars.gridVolVars().problem();
77 const Scalar forchCoeff = problem.spatialParams().forchCoeff(scvf);
87 static const Scalar epsilon = getParamFromGroup<Scalar>(problem.paramGroup(),
"Forchheimer.NewtonTolerance", 1e-12);
88 static const std::size_t maxNumIter = getParamFromGroup<std::size_t>(problem.paramGroup(),
"Forchheimer.MaxIterations", 30);
89 for (
int k = 0; residual.two_norm() > epsilon ; ++k)
92 DUNE_THROW(
NumericalProblem,
"could not determine forchheimer velocity within "
93 << k <<
" iterations");
96 forchheimerResidual_(residual,
106 forchheimerDerivative_(gradF,
115 gradF.solve(deltaV, residual);
121 auto upwindTerm = [phaseIdx](
const auto& volVars){
return volVars.mobility(phaseIdx); };
122 const Scalar forchheimerUpwindMobility = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm,
velocity * scvf.unitOuterNormal(), phaseIdx);
125 if (forchheimerUpwindMobility > 0.0)
126 velocity /= forchheimerUpwindMobility;
135 template<
class Problem,
class ElementVolumeVariables,
136 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
137 std::enable_if_t<scalarPerm, int> = 0>
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf)
143 Scalar harmonicMeanSqrtK(0.0);
145 const auto insideScvIdx = scvf.insideScvIdx();
146 const auto& insideVolVars = elemVolVars[insideScvIdx];
147 const Scalar Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
148 const Scalar sqrtKi = sqrt(Ki);
150 if (!scvf.boundary())
152 const auto outsideScvIdx = scvf.outsideScvIdx();
153 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
154 const Scalar Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
155 const Scalar sqrtKj = sqrt(Kj);
159 harmonicMeanSqrtK = sqrtKi;
163 for (
int i = 0; i < dimWorld; ++i)
164 result[i][i] = harmonicMeanSqrtK;
173 template<
class Problem,
class ElementVolumeVariables,
174 bool scalarPerm = std::is_same<typename Problem::SpatialParams::PermeabilityType, Scalar>::value,
175 std::enable_if_t<!scalarPerm, int> = 0>
177 const ElementVolumeVariables& elemVolVars,
178 const SubControlVolumeFace& scvf)
185 const auto insideScvIdx = scvf.insideScvIdx();
186 const auto& insideVolVars = elemVolVars[insideScvIdx];
187 const auto& Ki = getPermeability_(problem, insideVolVars, scvf.ipGlobal());
190 if (!isDiagonal_(Ki))
191 DUNE_THROW(Dune::NotImplemented,
"Only diagonal permeability tensors are supported.");
193 for (
int i = 0; i < dim; ++i)
194 sqrtKi[i][i] = sqrt(Ki[i][i]);
196 if (!scvf.boundary())
198 const auto outsideScvIdx = scvf.outsideScvIdx();
199 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
200 const auto& Kj = getPermeability_(problem, outsideVolVars, scvf.ipGlobal());
202 if (!isDiagonal_(Kj))
203 DUNE_THROW(Dune::NotImplemented,
"Only diagonal permeability tensors are supported.");
205 for (
int i = 0; i < dim; ++i)
206 sqrtKj[i][i] = sqrt(Kj[i][i]);
211 harmonicMeanSqrtK = sqrtKi;
213 return harmonicMeanSqrtK;
219 template<
class ElementVolumeVariables>
221 const Scalar forchCoeff,
225 const ElementVolumeVariables& elemVolVars,
226 const SubControlVolumeFace& scvf,
229 residual = oldVelocity;
231 residual -= darcyVelocity;
233 auto upwindTerm = [phaseIdx](
const auto& volVars){
return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
234 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm, oldVelocity * scvf.unitOuterNormal(), phaseIdx);
235 sqrtK.usmv(forchCoeff * rhoOverMu * oldVelocity.two_norm(), oldVelocity, residual);
239 template<
class ElementVolumeVariables>
241 const Scalar forchCoeff,
244 const ElementVolumeVariables& elemVolVars,
245 const SubControlVolumeFace& scvf,
251 auto upwindTerm = [phaseIdx](
const auto& volVars){
return volVars.density(phaseIdx) / volVars.viscosity(phaseIdx); };
253 const Scalar absV =
velocity.two_norm();
263 const Scalar rhoOverMu = UpwindScheme::multiplier(elemVolVars, scvf, upwindTerm,
velocity * scvf.unitOuterNormal(), phaseIdx);
266 for (
int i = 0; i < dim; i++)
268 for (
int k = 0; k <= i; k++)
274 derivative[k][i] = derivative[i][k];
281 for (
int i = 0; i < dim; i++)
282 derivative[i][i] += (1.0 + (sqrtK[i][i]*forchCoeff * absV * rhoOverMu));
296 for (
int i = 0; i < dim; i++)
297 for (
int k = 0; k < dim; k++)
298 if ((i != k) && (abs(K[i][k]) >= 1e-25))
303 template<
class Problem,
class VolumeVariables>
304 static decltype(
auto) getPermeability_(
const Problem& problem,
305 const VolumeVariables& volVars,
306 const GlobalPosition& scvfIpGlobal)
308 if constexpr (Problem::SpatialParams::evaluatePermeabilityAtScvfIP())
309 return problem.spatialParams().permeabilityAtPos(scvfIpGlobal);
311 return volVars.permeability();
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
Define some often used mathematical functions.
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
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
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 For a detailed description see dumux/flow/forchheimerslaw.hh.
Definition: forchheimervelocity.hh:51
typename FluxVariables::UpwindScheme UpwindScheme
Definition: forchheimervelocity.hh:61
static DimWorldMatrix calculateHarmonicMeanSqrtPermeability(const Problem &problem, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf)
Returns the harmonic mean of and .
Definition: forchheimervelocity.hh:138
Dune::FieldVector< Scalar, dimWorld > DimWorldVector
Definition: forchheimervelocity.hh:62
static DimWorldVector velocity(const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const DimWorldMatrix sqrtK, const DimWorldVector darcyVelocity)
Definition: forchheimervelocity.hh:67
Dune::FieldMatrix< Scalar, dimWorld, dimWorld > DimWorldMatrix
Definition: forchheimervelocity.hh:63
Declares all properties used in Dumux.