24#ifndef DUMUX_MATERIAL_FLUIDMATRIX_DISPERSIONTENSORS_SCHEIDEGGER_HH
25#define DUMUX_MATERIAL_FLUIDMATRIX_DISPERSIONTENSORS_SCHEIDEGGER_HH
29#include <dune/common/math.hh>
30#include <dune/common/std/type_traits.hh>
31#include <dune/common/fmatrix.hh>
41template <
class Problem,
class SubControlVolumeFace>
44template<
class Problem,
class SubControlVolumeFace>
46{
return Dune::Std::is_detected<HasVelocityInSpatialParams, Problem, SubControlVolumeFace>::value; }
53template<
class TypeTag>
58 using FVElementGeometry =
typename GridGeometry::LocalView;
59 using SubControlVolumeFace =
typename GridGeometry::SubControlVolumeFace;
64 static const int dimWorld = GridView::dimensionworld;
65 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
72 template <
class ElementFluxVariablesCache>
74 const SubControlVolumeFace& scvf,
75 const FVElementGeometry& fvGeometry,
76 const ElementVolumeVariables& elemVolVars,
77 const ElementFluxVariablesCache& elemFluxVarsCache,
81 DimWorldMatrix dispersionTensor(0.0);
84 auto velocity = dispersionVelocity_(problem, scvf, fvGeometry, elemVolVars, elemFluxVarsCache);
87 std::array<Scalar,2> dispersivity = problem.spatialParams().dispersionAlphas(scvf.center(), phaseIdx, compIdx);
89 return scheideggerTensor_(dispersivity, velocity);
92 template <
class ElementFluxVariablesCache>
94 const SubControlVolumeFace& scvf,
95 const FVElementGeometry& fvGeometry,
96 const ElementVolumeVariables& elemVolVars,
97 const ElementFluxVariablesCache& elemFluxVarsCache,
100 DimWorldMatrix dispersionTensor(0.0);
103 auto velocity = dispersionVelocity_(problem, scvf, fvGeometry, elemVolVars, elemFluxVarsCache);
106 std::array<Scalar,2> dispersivity = problem.spatialParams().dispersionAlphas(scvf.center(), phaseIdx);
108 return scheideggerTensor_(dispersivity, velocity);
113 template <
class ElementFluxVariablesCache>
114 static Dune::FieldVector<Scalar, dimWorld> dispersionVelocity_(
const Problem& problem,
115 const SubControlVolumeFace& scvf,
116 [[maybe_unused]]
const FVElementGeometry& fvGeometry,
117 [[maybe_unused]]
const ElementVolumeVariables& elemVolVars,
118 [[maybe_unused]]
const ElementFluxVariablesCache& elemFluxVarsCache)
121 Dune::FieldVector<Scalar, dimWorld> velocity(0.0);
122 if constexpr (stationaryVelocityField)
124 if constexpr (!Detail::hasVelocityInSpatialParams<Problem,SubControlVolumeFace>() )
125 DUNE_THROW(Dune::NotImplemented,
"\n Please provide the stationary velocity field in the spatialparams via a velocity function.");
127 velocity = problem.spatialParams().velocity(scvf);
133 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
134 const auto& shapeValues = fluxVarsCache.shapeValues();
137 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
138 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
140 outsideVolVars.permeability(),
141 scvf.unitOuterNormal());
144 Dune::FieldVector<Scalar, dimWorld> gradP(0.0);
146 static const bool enableGravity = getParamFromGroup<bool>(problem.paramGroup(),
"Problem.EnableGravity");
147 for (
auto&& scv : scvs(fvGeometry))
149 const auto& volVars = elemVolVars[scv];
152 rho += volVars.density(0)*shapeValues[scv.indexInElement()][0];
154 gradP.axpy(volVars.pressure(0), fluxVarsCache.gradN(scv.indexInElement()));
158 gradP.axpy(-rho, problem.spatialParams().gravity(scvf.center()));
163 velocity /= -0.5 * (insideVolVars.viscosity() + outsideVolVars.viscosity());
166 DUNE_THROW(Dune::NotImplemented,
"\n Scheidegger Dispersion for compositional models without given constant velocity field is only implemented using the Box method.");
172 static DimWorldMatrix scheideggerTensor_(
const std::array<Scalar,2>& dispersivity,
173 const Dune::FieldVector<Scalar, dimWorld>& velocity)
175 DimWorldMatrix scheideggerTensor(0.0);
178 for (
int i=0; i < dimWorld; i++)
179 for (
int j = 0; j < dimWorld; j++)
180 scheideggerTensor[i][j] = velocity[i]*velocity[j];
183 Scalar vNorm = velocity.two_norm();
185 scheideggerTensor /= vNorm;
187 scheideggerTensor = 0;
190 scheideggerTensor *= (dispersivity[0] - dispersivity[1]);
193 for (
int i = 0; i < dimWorld; i++)
194 scheideggerTensor[i][i] += vNorm*dispersivity[1];
196 return scheideggerTensor;
The infrastructure to retrieve run-time parameters from Dune::ParameterTrees.
The available discretization methods in Dumux.
A free function to average a Tensor at an interface.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:180
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
decltype(std::declval< Problem >().spatialParams().velocity(std::declval< SubControlVolumeFace >())) HasVelocityInSpatialParams
Definition: scheidegger.hh:42
static constexpr bool hasVelocityInSpatialParams()
Definition: scheidegger.hh:45
constexpr Box box
Definition: method.hh:136
Traits of a flux variables type.
Definition: flux/traits.hh:44
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:45
Scheidegger's dispersion tensor.
Definition: scheidegger.hh:55
static DimWorldMatrix compositionalDispersionTensor(const Problem &problem, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const int phaseIdx, const int compIdx)
Definition: scheidegger.hh:73
static DimWorldMatrix thermalDispersionTensor(const Problem &problem, const SubControlVolumeFace &scvf, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const ElementFluxVariablesCache &elemFluxVarsCache, const int phaseIdx)
Definition: scheidegger.hh:93
Declares all properties used in Dumux.