Line | Branch | Exec | Source |
---|---|---|---|
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 | // | ||
7 | /*! | ||
8 | * \file | ||
9 | * \ingroup CCTpfaFlux | ||
10 | * \brief This file contains the data which is required to calculate | ||
11 | * dispersive fluxes. | ||
12 | */ | ||
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> | ||
20 | #include <dumux/common/properties.hh> | ||
21 | #include <dumux/discretization/method.hh> | ||
22 | #include <dumux/discretization/extrusion.hh> | ||
23 | #include <dumux/discretization/cellcentered/tpfa/computetransmissibility.hh> | ||
24 | #include <dumux/flux/traits.hh> | ||
25 | #include <dumux/flux/referencesystemformulation.hh> | ||
26 | |||
27 | namespace Dumux { | ||
28 | |||
29 | // forward declaration | ||
30 | template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem> | ||
31 | class DispersionFluxImplementation; | ||
32 | |||
33 | /*! | ||
34 | * \ingroup CCTpfaFlux | ||
35 | * \brief Specialization of a Dispersion flux for the cctpfa method | ||
36 | */ | ||
37 | template <class TypeTag, ReferenceSystemFormulation referenceSystem> | ||
38 | class DispersionFluxImplementation<TypeTag, DiscretizationMethods::CCTpfa, referenceSystem> | ||
39 | { | ||
40 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
41 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
42 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
43 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
44 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
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; | ||
50 | using GridFluxVariablesCache = GetPropType<TypeTag, Properties::GridFluxVariablesCache>; | ||
51 | using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView; | ||
52 | using FluxVarCache = typename GridFluxVariablesCache::FluxVariablesCache; | ||
53 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
54 | using FluxTraits = typename Dumux::FluxTraits<FluxVariables>; | ||
55 | using BalanceEqOpts = GetPropType<TypeTag, Properties::BalanceEqOpts>; | ||
56 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
57 | using Element = typename GridView::template Codim<0>::Entity; | ||
58 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
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 | |||
73 | public: | ||
74 | |||
75 | //return the reference system | ||
76 | static constexpr ReferenceSystemFormulation referenceSystemFormulation() | ||
77 | { return referenceSystem; } | ||
78 | |||
79 | /*! | ||
80 | * \brief Returns the dispersive fluxes of all components within | ||
81 | * a fluid phase across the given sub-control volume face. | ||
82 | * The computed fluxes are given in mole/s or kg/s, depending | ||
83 | * on the template parameter ReferenceSystemFormulation. | ||
84 | */ | ||
85 | 989800 | 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 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 989800 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 989800 times.
|
1979600 | 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 | 989800 | ComponentFluxVector componentFlux(0.0); | |
99 | 1979600 | const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()]; | |
100 | 1979600 | const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()]; | |
101 | |||
102 | 989800 | const auto rhoInside = massOrMolarDensity(insideVolVars, referenceSystem, phaseIdx); | |
103 | 989800 | const auto rhoOutside = massOrMolarDensity(outsideVolVars, referenceSystem, phaseIdx); | |
104 | 989800 | const Scalar rho = 0.5*(rhoInside + rhoOutside); | |
105 | |||
106 |
2/2✓ Branch 0 taken 1979600 times.
✓ Branch 1 taken 989800 times.
|
2969400 | for (int compIdx = 0; compIdx < numComponents; compIdx++) |
107 | { | ||
108 | 1979600 | const auto& dispersionTensor = | |
109 | 1979600 | 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 | 1979600 | const auto tij = calculateTransmissibility_(fvGeometry, elemVolVars, scvf, dispersionTensor); | |
116 | |||
117 | 1979600 | const auto xInside = massOrMoleFraction(insideVolVars, referenceSystem, phaseIdx, compIdx); | |
118 | 1979600 | const auto xOutide = massOrMoleFraction(outsideVolVars, referenceSystem, phaseIdx, compIdx); | |
119 | |||
120 | 3959200 | componentFlux[compIdx] = (rho * tij * (xInside-xOutide)); | |
121 | } | ||
122 | 989800 | return componentFlux; | |
123 | } | ||
124 | |||
125 | /*! | ||
126 | * \brief Returns the thermal dispersive flux | ||
127 | * across the given sub-control volume face. | ||
128 | */ | ||
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 | |||
162 | private: | ||
163 | |||
164 | //! Compute transmissibilities | ||
165 | template<class Tensor> | ||
166 | 1979600 | static Scalar calculateTransmissibility_(const FVElementGeometry& fvGeometry, | |
167 | const ElementVolumeVariables& elemVolVars, | ||
168 | const SubControlVolumeFace& scvf, | ||
169 | const Tensor& D) | ||
170 | { | ||
171 | Scalar tij; | ||
172 | |||
173 | 1979600 | const auto insideScvIdx = scvf.insideScvIdx(); | |
174 | 1979600 | const auto& insideScv = fvGeometry.scv(insideScvIdx); | |
175 | 1979600 | const auto& insideVolVars = elemVolVars[insideScvIdx]; | |
176 | |||
177 | 1979600 | const Scalar ti = computeTpfaTransmissibility(fvGeometry, scvf, insideScv, D, insideVolVars.extrusionFactor()); | |
178 | |||
179 | // for the boundary (dirichlet) we only need ti | ||
180 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 1979600 times.
|
1979600 | if (scvf.boundary()) |
181 | { | ||
182 | ✗ | tij = Extrusion::area(fvGeometry, scvf)*ti; | |
183 | } | ||
184 | // otherwise we compute a tpfa harmonic mean | ||
185 | else | ||
186 | { | ||
187 | 1979600 | const auto outsideScvIdx = scvf.outsideScvIdx(); | |
188 | 1979600 | const auto& outsideScv = fvGeometry.scv(outsideScvIdx); | |
189 | 1979600 | 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 | 1979600 | 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 |
1/2✓ Branch 0 taken 1979600 times.
✗ Branch 1 not taken.
|
1979600 | if (ti*tj <= 0.0) |
200 | tij = 0; | ||
201 | else | ||
202 | 3959200 | tij = Extrusion::area(fvGeometry, scvf)*(ti * tj)/(ti + tj); | |
203 | } | ||
204 | |||
205 | 1979600 | return tij; | |
206 | } | ||
207 | |||
208 | }; | ||
209 | |||
210 | } // end namespace Dumux | ||
211 | |||
212 | #endif | ||
213 |