version 3.11-dev
All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Modules Pages Concepts
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-FileCopyrightText: 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 if constexpr (!FluidSystem::isTracerFluidSystem())
109 if (compIdx == FluidSystem::getMainComponent(phaseIdx))
110 continue;
111
112 const auto& dispersionTensor =
113 ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry,
114 elemVolVars, elemFluxVarsCache,
115 phaseIdx, compIdx);
116
117 // Here we assume that the Dispersion tensor is given at the scvfs such that Di = Dj
118 // This means that we don't have to distinguish between inside and outside tensors
119 const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor);
120
121 const auto xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx);
122 const auto xOutide = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx);
123
124 componentFlux[compIdx] = (rho * tij * (xInside-xOutide));
125 if constexpr (!FluidSystem::isTracerFluidSystem())
126 if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx))
127 componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx];
128 }
129 return componentFlux;
130 }
131
136 static HeatFluxScalar thermalDispersionFlux(const Problem& problem,
137 const Element& element,
138 const FVElementGeometry& fvGeometry,
139 const ElementVolumeVariables& elemVolVars,
140 const SubControlVolumeFace& scvf,
141 const int phaseIdx,
142 const ElementFluxVariablesCache& elemFluxVarsCache)
143 {
144 if (scvf.numOutsideScvs() > 1 )
145 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for conforming grids.");
146 if (!stationaryVelocityField)
147 DUNE_THROW(Dune::NotImplemented, "\n Dispersion using ccTPFA is only implemented for problems with stationary velocity fields");
148
149 const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()];
150 const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()];
151
152 const auto& dispersionTensor =
153 ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry,
154 elemVolVars, elemFluxVarsCache,
155 phaseIdx);
156
157 // Here we assume that the Dispersion tensor is given at the scvfs such that Di = Dj
158 // This means that we don't have to distinguish between inside and outside tensors
159 const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor);
160
161 // get the inside/outside temperatures
162 const auto tInside = insideVolVars.temperature();
163 const auto tOutside = outsideVolVars.temperature();
164
165 // compute the heat conduction flux
166 return tij * (tInside-tOutside);
167 }
168
169private:
170
172 template<class Tensor>
173 static Scalar calculateTransmissibility_(const FVElementGeometry& fvGeometry,
174 const ElementVolumeVariables& elemVolVars,
175 const SubControlVolumeFace& scvf,
176 const Tensor& D)
177 {
178 Scalar tij;
179
180 const auto insideScvIdx = scvf.insideScvIdx();
181 const auto& insideScv = fvGeometry.scv(insideScvIdx);
182 const auto& insideVolVars = elemVolVars[insideScvIdx];
183
184 const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv, D, insideVolVars.extrusionFactor());
185
186 // for the boundary (dirichlet) we only need ti
187 if (scvf.boundary())
188 {
189 tij = Extrusion::area(fvGeometry, scvf)*ti;
190 }
191 // otherwise we compute a tpfa harmonic mean
192 else
193 {
194 const auto outsideScvIdx = scvf.outsideScvIdx();
195 const auto& outsideScv = fvGeometry.scv(outsideScvIdx);
196 const auto& outsideVolVars = elemVolVars[outsideScvIdx];
197
198 Scalar tj;
199 if constexpr (dim == dimWorld)
200 // assume the normal vector from outside is anti parallel so we save flipping a vector
201 tj = -1.0*computeTpfaTransmissibility(fvGeometry, scvf, outsideScv, D, outsideVolVars.extrusionFactor());
202 else
203 tj = computeTpfaTransmissibility(fvGeometry, fvGeometry.flipScvf(scvf.index()), outsideScv, D, outsideVolVars.extrusionFactor());
204
205 // check for division by zero!
206 if (ti*tj <= 0.0)
207 tij = 0;
208 else
209 tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj);
210 }
211
212 return tij;
213 }
214
215};
216
217} // end namespace Dumux
218
219#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:136
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...