version 3.10-dev
box/dispersionflux.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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
13#ifndef DUMUX_DISCRETIZATION_BOX_DISPERSION_FLUX_HH
14#define DUMUX_DISCRETIZATION_BOX_DISPERSION_FLUX_HH
15
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
18#include <dune/common/float_cmp.hh>
19
20#include <dumux/common/math.hh>
24#include <dumux/flux/traits.hh>
27
28namespace Dumux {
29
30// forward declaration
31template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
33
38template <class TypeTag, ReferenceSystemFormulation referenceSystem>
39class DispersionFluxImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem>
40{
46 using FVElementGeometry = typename GridGeometry::LocalView;
47 using SubControlVolume = typename GridGeometry::SubControlVolume;
48 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
49 using Extrusion = Extrusion_t<GridGeometry>;
50 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
52 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
53 using FluxVarCache = typename GridFluxVariablesCache::FluxVariablesCache;
55 using FluxTraits = typename Dumux::FluxTraits<FluxVariables>;
58 using Element = typename GridView::template Codim<0>::Entity;
60 using Indices = typename ModelTraits::Indices;
61
62 static constexpr int dim = GridView::dimension;
63 static constexpr int dimWorld = GridView::dimensionworld;
64
65 static constexpr int numPhases = ModelTraits::numFluidPhases();
66 static constexpr int numComponents = ModelTraits::numFluidComponents();
67
68 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
69 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
70 using HeatFluxScalar = Scalar;
71
72 static constexpr bool stationaryVelocityField = FluxTraits::hasStationaryVelocityField();
73
74public:
75
76 //return the reference system
78 { return referenceSystem; }
79
86 static ComponentFluxVector compositionalDispersionFlux(const Problem& problem,
87 const Element& element,
88 const FVElementGeometry& fvGeometry,
89 const ElementVolumeVariables& elemVolVars,
90 const SubControlVolumeFace& scvf,
91 const int phaseIdx,
92 const ElementFluxVariablesCache& elemFluxVarsCache)
93 {
94 ComponentFluxVector componentFlux(0.0);
95
96 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
97 const auto& shapeValues = fluxVarsCache.shapeValues();
98
99 // density interpolation
100 Scalar rhoMassOrMole(0.0);
101 for (auto&& scv : scvs(fvGeometry))
102 {
103 const auto rho = massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx);
104 rhoMassOrMole += rho * shapeValues[scv.indexInElement()][0];
105 }
106
107 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
108 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
109
110 for (int compIdx = 0; compIdx < numComponents; compIdx++)
111 {
112 // collect the dispersion tensor, the fluxVarsCache and the shape values
113 const auto& dispersionTensor = [&]()
114 {
115 const auto& tensor =
116 ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
117 elemVolVars, elemFluxVarsCache,
118 phaseIdx, compIdx);
119
120 if (Dune::FloatCmp::eq(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor(), 1e-6))
121 return insideVolVars.extrusionFactor()*tensor;
122 else
123 {
124 const auto insideTensor = insideVolVars.extrusionFactor() * tensor;
125 const auto outsideTensor = outsideVolVars.extrusionFactor() * tensor;
126
127 // the resulting averaged dispersion tensor
128 return faceTensorAverage(insideTensor, outsideTensor, scvf.unitOuterNormal());
129 }
130 }();
131
132 // the mole/mass fraction gradient
133 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
134 for (auto&& scv : scvs(fvGeometry))
135 {
136 const auto x = massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx);
137 gradX.axpy(x, fluxVarsCache.gradN(scv.indexInElement()));
138 }
139
140 // compute the dispersion flux
141 componentFlux[compIdx] = -1.0 * rhoMassOrMole * vtmv(scvf.unitOuterNormal(), dispersionTensor, gradX)*Extrusion::area(fvGeometry, scvf);
142 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
143 componentFlux[phaseIdx] -= componentFlux[compIdx];
144 }
145 return componentFlux;
146 }
147
152 static HeatFluxScalar thermalDispersionFlux(const Problem& problem,
153 const Element& element,
154 const FVElementGeometry& fvGeometry,
155 const ElementVolumeVariables& elemVolVars,
156 const SubControlVolumeFace& scvf,
157 const int phaseIdx,
158 const ElementFluxVariablesCache& elemFluxVarsCache)
159 {
160 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
161 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
162
163 // collect the dispersion tensor
164 const auto& dispersionTensor = [&]()
165 {
166 const auto& tensor =
167 ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
168 elemVolVars, elemFluxVarsCache,
169 phaseIdx);
170
171 if (Dune::FloatCmp::eq(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor(), 1e-6))
172 return insideVolVars.extrusionFactor()*tensor;
173 else
174 {
175 const auto insideTensor = insideVolVars.extrusionFactor() * tensor;
176 const auto outsideTensor = outsideVolVars.extrusionFactor() * tensor;
177
178 // the resulting averaged dispersion tensor
179 return faceTensorAverage(insideTensor, outsideTensor, scvf.unitOuterNormal());
180 }
181 }();
182
183 // compute the temperature gradient with the shape functions
184 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
185 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
186 for (auto&& scv : scvs(fvGeometry))
187 gradTemp.axpy(elemVolVars[scv].temperature(), fluxVarsCache.gradN(scv.indexInElement()));
188
189 // compute the heat conduction flux
190 return -1.0*vtmv(scvf.unitOuterNormal(), dispersionTensor, gradTemp)*Extrusion::area(fvGeometry, scvf);
191 }
192
193};
194
195} // end namespace Dumux
196
197#endif
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/dispersionflux.hh:77
static ComponentFluxVector compositionalDispersionFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the dispersive fluxes of all components within a fluid phase across the given sub-control vol...
Definition: box/dispersionflux.hh:86
static HeatFluxScalar thermalDispersionFlux(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry, const ElementVolumeVariables &elemVolVars, const SubControlVolumeFace &scvf, const int phaseIdx, const ElementFluxVariablesCache &elemFluxVarsCache)
Returns the thermal dispersive flux across the given sub-control volume face.
Definition: box/dispersionflux.hh:152
Definition: box/dispersionflux.hh:32
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
A free function to average a Tensor at an interface.
Defines the flux traits.
Dune::DenseMatrix< MAT >::value_type vtmv(const Dune::DenseVector< V1 > &v1, const Dune::DenseMatrix< MAT > &M, const Dune::DenseVector< V2 > &v2)
Evaluates the scalar product of a vector v2, projected by a matrix M, with a vector v1.
Definition: math.hh:880
VolumeVariables::PrimaryVariables::value_type massOrMoleFraction(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx, const int compIdx)
returns the mass or mole fraction to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:54
VolumeVariables::PrimaryVariables::value_type massOrMolarDensity(const VolumeVariables &volVars, ReferenceSystemFormulation referenceSys, const int phaseIdx)
evaluates the density to be used in Fick's law based on the reference system
Definition: referencesystemformulation.hh:43
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:33
typename GetProp< TypeTag, Property >::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:296
Define some often used mathematical functions.
The available discretization methods in Dumux.
CVFE< CVFEMethods::PQ1 > Box
Definition: method.hh:94
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:39
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
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:29
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
Traits of a flux variables type.
Definition: flux/traits.hh:32
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:33