version 3.9
cctpfa/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_CC_TPFA_DISPERSION_FLUX_HH
14#define DUMUX_DISCRETIZATION_CC_TPFA_DISPERSION_FLUX_HH
15
16#include <dune/common/fvector.hh>
17#include <dune/common/fmatrix.hh>
18
19#include <dumux/common/math.hh>
24#include <dumux/flux/traits.hh>
26
27namespace Dumux {
28
29// forward declaration
30template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem>
31class DispersionFluxImplementation;
32
37template <class TypeTag, ReferenceSystemFormulation referenceSystem>
38class DispersionFluxImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem>
39{
45 using FVElementGeometry = typename GridGeometry::LocalView;
46 using SubControlVolume = typename GridGeometry::SubControlVolume;
47 using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace;
48 using Extrusion = Extrusion_t<GridGeometry>;
49 using ElementVolumeVariables = typename GetPropType<TypeTag, Properties::GridVolumeVariables>::LocalView;
51 using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView;
52 using FluxVarCache = typename GridFluxVariablesCache::FluxVariablesCache;
54 using FluxTraits = typename Dumux::FluxTraits<FluxVariables>;
57 using Element = typename GridView::template Codim<0>::Entity;
59 using Indices = typename ModelTraits::Indices;
60
61 static constexpr int dim = GridView::dimension;
62 static constexpr int dimWorld = GridView::dimensionworld;
63
64 static constexpr int numPhases = ModelTraits::numFluidPhases();
65 static constexpr int numComponents = ModelTraits::numFluidComponents();
66
67 using DimWorldMatrix = Dune::FieldMatrix<Scalar, dimWorld, dimWorld>;
68 using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>;
69 using HeatFluxScalar = Scalar;
70
71 static constexpr bool stationaryVelocityField = FluxTraits::hasStationaryVelocityField();
72
73public:
74
75 //return the reference system
77 { return referenceSystem; }
78
85 static ComponentFluxVector compositionalDispersionFlux(const Problem& problem,
86 const Element& element,
87 const FVElementGeometry& fvGeometry,
88 const ElementVolumeVariables& elemVolVars,
89 const SubControlVolumeFace& scvf,
90 const int phaseIdx,
91 const ElementFluxVariablesCache& elemFluxVarsCache)
92 {
93 if (scvf.numOutsideScvs() > 1 )
94 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for conforming grids.");
95 if (!stationaryVelocityField)
96 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
97
98 ComponentFluxVector componentFlux(0.0);
99 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
100 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
101
102 const auto rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx);
103 const auto rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx);
104 const Scalar rho = 0.5*(rhoInside + rhoOutside);
105
106 for (int compIdx = 0; compIdx < numComponents; compIdx++)
107 {
108 const auto& dispersionTensor =
109 ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
110 elemVolVars, elemFluxVarsCache,
111 phaseIdx, compIdx);
112
113 // Here we assume that the Dispersion tensor is given at the scvfs such that Di = Dj
114 // This means that we don't have to distinguish between inside and outside tensors
115 const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor);
116
117 const auto xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
118 const auto xOutide = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
119
120 componentFlux[compIdx] = (rho * tij * (xInside-xOutide));
121 }
122 return componentFlux;
123 }
124
129 static HeatFluxScalar thermalDispersionFlux(const Problem& problem,
130 const Element& element,
131 const FVElementGeometry& fvGeometry,
132 const ElementVolumeVariables& elemVolVars,
133 const SubControlVolumeFace& scvf,
134 const int phaseIdx,
135 const ElementFluxVariablesCache& elemFluxVarsCache)
136 {
137 if (scvf.numOutsideScvs() > 1 )
138 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for conforming grids.");
139 if (!stationaryVelocityField)
140 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
141
142 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
143 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
144
145 const auto& dispersionTensor =
146 ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
147 elemVolVars, elemFluxVarsCache,
148 phaseIdx);
149
150 // Here we assume that the Dispersion tensor is given at the scvfs such that Di = Dj
151 // This means that we don't have to distinguish between inside and outside tensors
152 const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor);
153
154 // get the inside/outside temperatures
155 const auto tInside = insideVolVars.temperature();
156 const auto tOutside = outsideVolVars.temperature();
157
158 // compute the heat conduction flux
159 return tij * (tInside-tOutside);
160 }
161
162private:
163
165 template<class Tensor>
166 static Scalar calculateTransmissibility_(const FVElementGeometry& fvGeometry,
167 const ElementVolumeVariables& elemVolVars,
168 const SubControlVolumeFace& scvf,
169 const Tensor& D)
170 {
171 Scalar tij;
172
173 const auto insideScvIdx = scvf.insideScvIdx();
174 const auto& insideScv = fvGeometry.scv(insideScvIdx);
175 const auto& insideVolVars = elemVolVars[insideScvIdx];
176
177 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv, D, insideVolVars.extrusionFactor());
178
179 // for the boundary (dirichlet) we only need ti
180 if (scvf.boundary())
181 {
182 tij = Extrusion::area(fvGeometry, scvf)*ti;
183 }
184 // otherwise we compute a tpfa harmonic mean
185 else
186 {
187 const auto outsideScvIdx = scvf.outsideScvIdx();
188 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
189 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
190
191 Scalar tj;
192 if constexpr (dim == dimWorld)
193 // assume the normal vector from outside is anti parallel so we save flipping a vector
194 tj = -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, D, outsideVolVars.extrusionFactor());
195 else
196 tj = computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv, D, outsideVolVars.extrusionFactor());
197
198 // check for division by zero!
199 if (ti*tj <= 0.0)
200 tij = 0;
201 else
202 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
203 }
204
205 return tij;
206 }
207
208};
209
210} // end namespace Dumux
211
212#endif
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: cctpfa/dispersionflux.hh:85
static constexpr ReferenceSystemFormulation referenceSystemFormulation()
Definition: cctpfa/dispersionflux.hh:76
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: cctpfa/dispersionflux.hh:129
Definition: box/dispersionflux.hh:32
Defines all properties used in Dumux.
Helper classes to compute the integration elements.
Defines the flux traits.
Tensor::field_type computeTpfaTransmissibility(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf, const typename FVElementGeometry::SubControlVolume &scv, const Tensor &T, typename FVElementGeometry::SubControlVolume::Traits::Scalar extrusionFactor)
Free function to evaluate the Tpfa transmissibility associated with the flux (in the form of flux = T...
Definition: tpfa/computetransmissibility.hh:36
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.
Definition: adapt.hh:17
typename Extrusion< T >::type Extrusion_t
Convenience alias for obtaining the extrusion type.
Definition: extrusion.hh:166
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
Free functions to evaluate the transmissibilities associated with flux evaluations across sub-control...