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 BoxFlux | ||
10 | * \brief This file contains the data which is required to calculate | ||
11 | * dispersive fluxes. | ||
12 | */ | ||
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> | ||
21 | #include <dumux/common/properties.hh> | ||
22 | #include <dumux/discretization/method.hh> | ||
23 | #include <dumux/discretization/extrusion.hh> | ||
24 | #include <dumux/flux/traits.hh> | ||
25 | #include <dumux/flux/referencesystemformulation.hh> | ||
26 | #include <dumux/flux/facetensoraverage.hh> | ||
27 | |||
28 | namespace Dumux { | ||
29 | |||
30 | // forward declaration | ||
31 | template<class TypeTag, class DiscretizationMethod, ReferenceSystemFormulation referenceSystem> | ||
32 | class DispersionFluxImplementation; | ||
33 | |||
34 | /*! | ||
35 | * \ingroup BoxFlux | ||
36 | * \brief Specialization of a dispersion flux for the box method | ||
37 | */ | ||
38 | template <class TypeTag, ReferenceSystemFormulation referenceSystem> | ||
39 | class DispersionFluxImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem> | ||
40 | { | ||
41 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
42 | using Problem = GetPropType<TypeTag, Properties::Problem>; | ||
43 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
44 | using VolumeVariables = GetPropType<TypeTag, Properties::VolumeVariables>; | ||
45 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
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; | ||
51 | using GridFluxVariablesCache = GetPropType<TypeTag, Properties::GridFluxVariablesCache>; | ||
52 | using ElementFluxVariablesCache = typename GridFluxVariablesCache::LocalView; | ||
53 | using FluxVarCache = typename GridFluxVariablesCache::FluxVariablesCache; | ||
54 | using FluxVariables = GetPropType<TypeTag, Properties::FluxVariables>; | ||
55 | using FluxTraits = typename Dumux::FluxTraits<FluxVariables>; | ||
56 | using BalanceEqOpts = GetPropType<TypeTag, Properties::BalanceEqOpts>; | ||
57 | using GridView = typename GetPropType<TypeTag, Properties::GridGeometry>::GridView; | ||
58 | using Element = typename GridView::template Codim<0>::Entity; | ||
59 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
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 | |||
74 | public: | ||
75 | |||
76 | //return the reference system | ||
77 | static constexpr ReferenceSystemFormulation referenceSystemFormulation() | ||
78 | { return referenceSystem; } | ||
79 | |||
80 | /*! | ||
81 | * \brief Returns the dispersive fluxes of all components within | ||
82 | * a fluid phase across the given sub-control volume face. | ||
83 | * The computed fluxes are given in mole/s or kg/s, depending | ||
84 | * on the template parameter ReferenceSystemFormulation. | ||
85 | */ | ||
86 | 9667600 | 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 | 9667600 | ComponentFluxVector componentFlux(0.0); | |
95 | |||
96 | 9667600 | const auto& fluxVarsCache = elemFluxVarsCache[scvf]; | |
97 | 9667600 | const auto& shapeValues = fluxVarsCache.shapeValues(); | |
98 | |||
99 | // density interpolation | ||
100 | 9667600 | Scalar rhoMassOrMole(0.0); | |
101 |
4/4✓ Branch 0 taken 38670400 times.
✓ Branch 1 taken 9667600 times.
✓ Branch 2 taken 38670400 times.
✓ Branch 3 taken 9667600 times.
|
96676000 | for (auto&& scv : scvs(fvGeometry)) |
102 | { | ||
103 | 77340800 | const auto rho = massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx); | |
104 | 77340800 | rhoMassOrMole += rho * shapeValues[scv.indexInElement()][0]; | |
105 | } | ||
106 | |||
107 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 9667600 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 9667600 times.
|
19335200 | const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()]; |
108 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 9667600 times.
|
9667600 | const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()]; |
109 | |||
110 |
2/2✓ Branch 0 taken 19335200 times.
✓ Branch 1 taken 9667600 times.
|
29002800 | for (int compIdx = 0; compIdx < numComponents; compIdx++) |
111 | { | ||
112 | // collect the dispersion tensor, the fluxVarsCache and the shape values | ||
113 | 38670400 | const auto& dispersionTensor = [&]() | |
114 | { | ||
115 | 19335200 | const auto& tensor = | |
116 | 69797600 | ModelTraits::CompositionalDispersionModel::compositionalDispersionTensor(problem, scvf, fvGeometry, | |
117 | 19335200 | elemVolVars, elemFluxVarsCache, | |
118 | phaseIdx, compIdx); | ||
119 | |||
120 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 19335200 times.
✓ Branch 2 taken 19335200 times.
✗ Branch 3 not taken.
|
19335200 | if (Dune::FloatCmp::eq(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor(), 1e-6)) |
121 | 19335200 | 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 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 11792000 times.
|
38670400 | }(); |
131 | |||
132 | // the mole/mass fraction gradient | ||
133 | 19335200 | Dune::FieldVector<Scalar, dimWorld> gradX(0.0); | |
134 |
4/4✓ Branch 0 taken 77340800 times.
✓ Branch 1 taken 19335200 times.
✓ Branch 2 taken 77340800 times.
✓ Branch 3 taken 19335200 times.
|
193352000 | for (auto&& scv : scvs(fvGeometry)) |
135 | { | ||
136 | 154681600 | const auto x = massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); | |
137 | 232022400 | gradX.axpy(x, fluxVarsCache.gradN(scv.indexInElement())); | |
138 | } | ||
139 | |||
140 | // compute the dispersion flux | ||
141 | 77340800 | componentFlux[compIdx] = -1.0 * rhoMassOrMole * vtmv(scvf.unitOuterNormal(), dispersionTensor, gradX)*Extrusion::area(fvGeometry, scvf); | |
142 | 19335200 | if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx) && !FluidSystem::isTracerFluidSystem()) | |
143 | 51945600 | componentFlux[phaseIdx] -= componentFlux[compIdx]; | |
144 | } | ||
145 | 9667600 | return componentFlux; | |
146 | } | ||
147 | |||
148 | /*! | ||
149 | * \brief Returns the thermal dispersive flux | ||
150 | * across the given sub-control volume face. | ||
151 | */ | ||
152 | 5600000 | 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 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 5600000 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 5600000 times.
|
11200000 | const auto& insideVolVars = elemVolVars[scvf.insideScvIdx()]; |
161 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 5600000 times.
|
5600000 | const auto& outsideVolVars = elemVolVars[scvf.outsideScvIdx()]; |
162 | |||
163 | // collect the dispersion tensor | ||
164 | 11200000 | const auto& dispersionTensor = [&]() | |
165 | { | ||
166 | 5600000 | const auto& tensor = | |
167 | 14918400 | ModelTraits::ThermalDispersionModel::thermalDispersionTensor(problem, scvf, fvGeometry, | |
168 | 5600000 | elemVolVars, elemFluxVarsCache, | |
169 | phaseIdx); | ||
170 | |||
171 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 5600000 times.
✓ Branch 2 taken 5600000 times.
✗ Branch 3 not taken.
|
5600000 | if (Dune::FloatCmp::eq(insideVolVars.extrusionFactor(), outsideVolVars.extrusionFactor(), 1e-6)) |
172 | 5600000 | 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 |
1/2✗ Branch 1 not taken.
✓ Branch 2 taken 3718400 times.
|
11200000 | }(); |
182 | |||
183 | // compute the temperature gradient with the shape functions | ||
184 | 5600000 | const auto& fluxVarsCache = elemFluxVarsCache[scvf]; | |
185 | 5600000 | Dune::FieldVector<Scalar, GridView::dimensionworld> gradTemp(0.0); | |
186 |
4/4✓ Branch 0 taken 22400000 times.
✓ Branch 1 taken 5600000 times.
✓ Branch 2 taken 22400000 times.
✓ Branch 3 taken 5600000 times.
|
56000000 | for (auto&& scv : scvs(fvGeometry)) |
187 | 112000000 | gradTemp.axpy(elemVolVars[scv].temperature(), fluxVarsCache.gradN(scv.indexInElement())); | |
188 | |||
189 | // compute the heat conduction flux | ||
190 | 11200000 | return -1.0*vtmv(scvf.unitOuterNormal(), dispersionTensor, gradTemp)*Extrusion::area(fvGeometry, scvf); | |
191 | } | ||
192 | |||
193 | }; | ||
194 | |||
195 | } // end namespace Dumux | ||
196 | |||
197 | #endif | ||
198 |