3.5-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
25#ifndef DUMUX_DISCRETIZATION_BOX_DISPERSION_FLUX_HH
26#define DUMUX_DISCRETIZATION_BOX_DISPERSION_FLUX_HH
27
28#include <dune/common/fvector.hh>
29#include <dune/common/fmatrix.hh>
30
31#include <dumux/common/math.hh>
35#include <dumux/flux/traits.hh>
37
38namespace Dumux {
39
40// forward declaration
41template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
43
48template <class TypeTag, ReferenceSystemFormulation referenceSystem>
49class DispersionFluxImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem>
50{
56 using FVElementGeometry = typename GridGeometry::LocalView;
57 using SubControlVolume = typename GridGeometry::SubControlVolume;
58 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
59 using Extrusion = Extrusion_t<GridGeometry>;
60 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
61 using ElementFluxVariablesCache = typename GetPropType<TypeTag, Properties::GridFluxVariablesCache>::LocalView;
64 using FluxTraits = typename Dumux::FluxTraits<FluxVariables>;
67 using Element = typename GridView::template Codim<0>::Entity;
69 using Indices = typename ModelTraits::Indices;
70
71 enum { dim = GridView::dimension} ;
72 enum { dimWorld = GridView::dimensionworld} ;
73 enum
74 {
75 numPhases = ModelTraits::numFluidPhases(),
76 numComponents = ModelTraits::numFluidComponents()
77 };
78
79 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
80 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
81 using HeatFluxScalar = Scalar;
82
83 static constexpr bool stationaryVelocityField = FluxTraits::hasStationaryVelocityField();
84
85public:
86
87 //return the reference system
89 { return referenceSystem; }
90
97 static ComponentFluxVector compositionalDispersionFlux(const Problem& problem,
98 const Element& element,
99 const FVElementGeometry& fvGeometry,
100 const ElementVolumeVariables& elemVolVars,
101 const SubControlVolumeFace& scvf,
102 const int phaseIdx,
103 const ElementFluxVariablesCache& elemFluxVarsCache)
104 {
105 ComponentFluxVector componentFlux(0.0);
106
107 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
108 const auto& shapeValues = fluxVarsCache.shapeValues();
109
110 // density interpolation
111 Scalar rhoMassOrMole(0.0);
112 for (auto&& scv : scvs(fvGeometry))
113 {
114 const auto rho = massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx);
115 rhoMassOrMole += rho * shapeValues[scv.indexInElement()][0];
116 }
117
118 for (int compIdx = 0; compIdx < numComponents; compIdx++)
119 {
120 // collect the dispersion tensor, the fluxVarsCache and the shape values
121 const auto& dispersionTensor =
122 ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
123 elemVolVars, elemFluxVarsCache,
124 phaseIdx, compIdx);
125
126 // the mole/mass fraction gradient
127 Dune::FieldVector<Scalar, dimWorld> gradX(0.0);
128 for (auto&& scv : scvs(fvGeometry))
129 {
130 const auto x = massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx);
131 gradX.axpy(x, fluxVarsCache.gradN(scv.indexInElement()));
132 }
133
134 // compute the dispersion flux
135 componentFlux[compIdx] = -1.0 * rhoMassOrMole * vtmv(scvf.unitOuterNormal(), dispersionTensor, gradX) * scvf.area();
136 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem())
137 componentFlux[phaseIdx] -= componentFlux[compIdx];
138 }
139 return componentFlux;
140 }
141
146 static HeatFluxScalar thermalDispersionFlux(const Problem& problem,
147 const Element& element,
148 const FVElementGeometry& fvGeometry,
149 const ElementVolumeVariables& elemVolVars,
150 const SubControlVolumeFace& scvf,
151 const int phaseIdx,
152 const ElementFluxVariablesCache& elemFluxVarsCache)
153 {
154 // collect the dispersion tensor
155 const auto& dispersionTensor =
156 ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
157 elemVolVars, elemFluxVarsCache,
158 phaseIdx);
159 // compute the temperature gradient with the shape functions
160 const auto& fluxVarsCache = elemFluxVarsCache[scvf];
161 Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0);
162 for (auto&& scv : scvs(fvGeometry))
163 gradTemp.axpy(elemVolVars[scv].temperature(), fluxVarsCache.gradN(scv.indexInElement()));
164
165 // compute the heat conduction flux
166 return -1.0*vtmv(scvf.unitOuterNormal(), dispersionTensor, gradTemp)*Extrusion::area(scvf);
167 }
168
169};
170
171} // end namespace Dumux
172
173#endif
Define some often used mathematical functions.
Helper classes to compute the integration elements.
The available discretization methods in Dumux.
The reference frameworks and formulations available for splitting total fluxes into a advective and d...
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:66
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:55
ReferenceSystemFormulation
The formulations available for Fick's law related to the reference system.
Definition: referencesystemformulation.hh:45
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:863
Definition: adapt.hh:29
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:177
typename Properties::Detail::GetPropImpl< TypeTag, Property >::type::type GetPropType
get the type alias defined in the property
Definition: propertysystem.hh:150
std::string temperature() noexcept
I/O name of temperature for equilibrium models.
Definition: name.hh:51
Definition: box/dispersionflux.hh:42
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: box/dispersionflux.hh:88
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:97
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:146
Traits of a flux variables type.
Definition: flux/traits.hh:44
static constexpr bool hasStationaryVelocityField()
Definition: flux/traits.hh:45
Declares all properties used in Dumux.
Defines the flux traits.