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 FacetCoupling | ||
10 | * \copydoc Dumux::BoxFacetCouplingFicksLaw | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FICKS_LAW_HH | ||
13 | #define DUMUX_DISCRETIZATION_BOX_FACET_COUPLING_FICKS_LAW_HH | ||
14 | |||
15 | #include <vector> | ||
16 | #include <cmath> | ||
17 | |||
18 | #include <dune/common/exceptions.hh> | ||
19 | #include <dune/common/fvector.hh> | ||
20 | #include <dune/common/float_cmp.hh> | ||
21 | |||
22 | #include <dumux/common/parameters.hh> | ||
23 | #include <dumux/common/properties.hh> | ||
24 | |||
25 | #include <dumux/flux/referencesystemformulation.hh> | ||
26 | #include <dumux/flux/box/fickslaw.hh> | ||
27 | #include <dumux/discretization/method.hh> | ||
28 | #include <dumux/discretization/extrusion.hh> | ||
29 | |||
30 | namespace Dumux { | ||
31 | |||
32 | /*! | ||
33 | * \ingroup FacetCoupling | ||
34 | * \brief Ficks's law for the box scheme in the context of coupled models | ||
35 | * where coupling occurs across the facets of the bulk domain elements | ||
36 | * with a lower-dimensional domain living on these facets. | ||
37 | */ | ||
38 | template<class TypeTag, ReferenceSystemFormulation referenceSystem = ReferenceSystemFormulation::massAveraged> | ||
39 | class BoxFacetCouplingFicksLaw | ||
40 | : public FicksLawImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem> | ||
41 | { | ||
42 | using ParentType = FicksLawImplementation<TypeTag, DiscretizationMethods::Box, referenceSystem>; | ||
43 | |||
44 | using GridGeometry = GetPropType<TypeTag, Properties::GridGeometry>; | ||
45 | using FVElementGeometry = typename GridGeometry::LocalView; | ||
46 | using SubControlVolumeFace = typename GridGeometry::SubControlVolumeFace; | ||
47 | using Extrusion = Extrusion_t<GridGeometry>; | ||
48 | using GridView = typename GridGeometry::GridView; | ||
49 | using Element = typename GridView::template Codim<0>::Entity; | ||
50 | |||
51 | static constexpr int dim = GridView::dimension; | ||
52 | static constexpr int dimWorld = GridView::dimensionworld; | ||
53 | |||
54 | using FluidSystem = GetPropType<TypeTag, Properties::FluidSystem>; | ||
55 | using ModelTraits = GetPropType<TypeTag, Properties::ModelTraits>; | ||
56 | using BalanceEqOpts = GetPropType<TypeTag, Properties::BalanceEqOpts>; | ||
57 | static constexpr int numComponents = ModelTraits::numFluidComponents(); | ||
58 | |||
59 | using Scalar = GetPropType<TypeTag, Properties::Scalar>; | ||
60 | using ComponentFluxVector = Dune::FieldVector<Scalar, numComponents>; | ||
61 | |||
62 | public: | ||
63 | |||
64 | template<class Problem, class ElementVolumeVariables, class ElementFluxVarsCache> | ||
65 | 7849832 | static ComponentFluxVector flux(const Problem& problem, | |
66 | const Element& element, | ||
67 | const FVElementGeometry& fvGeometry, | ||
68 | const ElementVolumeVariables& elemVolVars, | ||
69 | const SubControlVolumeFace& scvf, | ||
70 | const int phaseIdx, | ||
71 | const ElementFluxVarsCache& elemFluxVarCache) | ||
72 | { | ||
73 | // if this scvf is not on an interior boundary, use the standard law | ||
74 |
2/2✓ Branch 0 taken 7201376 times.
✓ Branch 1 taken 648456 times.
|
7849832 | if (!scvf.interiorBoundary()) |
75 | 7201376 | return ParentType::flux(problem, element, fvGeometry, elemVolVars, scvf, phaseIdx, elemFluxVarCache); | |
76 | |||
77 |
5/8✓ Branch 0 taken 5 times.
✓ Branch 1 taken 648451 times.
✓ Branch 3 taken 5 times.
✗ Branch 4 not taken.
✓ Branch 6 taken 5 times.
✗ Branch 7 not taken.
✓ Branch 9 taken 5 times.
✗ Branch 10 not taken.
|
648456 | static const Scalar xi = getParamFromGroup<Scalar>(problem.paramGroup(), "FacetCoupling.Xi", 1.0); |
78 |
2/4✗ Branch 0 not taken.
✓ Branch 1 taken 648456 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 648456 times.
|
648456 | if ( !Dune::FloatCmp::eq(xi, 1.0, 1e-6) ) |
79 | ✗ | DUNE_THROW(Dune::NotImplemented, "Xi != 1.0 cannot be used with the Box-Facet-Coupling scheme"); | |
80 | |||
81 | // get some references for convenience | ||
82 | 648456 | const auto& fluxVarCache = elemFluxVarCache[scvf]; | |
83 | 648456 | const auto& shapeValues = fluxVarCache.shapeValues(); | |
84 | 1296912 | const auto& insideScv = fvGeometry.scv(scvf.insideScvIdx()); | |
85 | 648456 | const auto& insideVolVars = elemVolVars[insideScv]; | |
86 | |||
87 | // interpolate density to scvf integration point | ||
88 | 648456 | Scalar rho = 0.0; | |
89 |
4/4✓ Branch 0 taken 2593824 times.
✓ Branch 1 taken 648456 times.
✓ Branch 2 taken 2593824 times.
✓ Branch 3 taken 648456 times.
|
6484560 | for (const auto& scv : scvs(fvGeometry)) |
90 | 10375296 | rho += massOrMolarDensity(elemVolVars[scv], referenceSystem, phaseIdx)*shapeValues[scv.indexInElement()][0]; | |
91 | |||
92 | // on interior Neumann boundaries, evaluate the flux using the facet effective diffusion coefficient | ||
93 | 648456 | ComponentFluxVector componentFlux(0.0); | |
94 | 648456 | const auto bcTypes = problem.interiorBoundaryTypes(element, scvf); | |
95 |
2/4✓ Branch 0 taken 648456 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 648456 times.
✗ Branch 3 not taken.
|
1296912 | if (bcTypes.hasOnlyNeumann()) |
96 | { | ||
97 | // compute tpfa flux from integration point to facet centerline | ||
98 | 1296912 | const auto& facetVolVars = problem.couplingManager().getLowDimVolVars(element, scvf); | |
99 | |||
100 | using std::sqrt; | ||
101 | // If this is a surface grid, use the square root of the facet extrusion factor | ||
102 | // as an approximate average distance from scvf ip to facet center | ||
103 | using std::sqrt; | ||
104 | 648456 | const auto a = facetVolVars.extrusionFactor(); | |
105 | 648456 | auto preGradX = scvf.unitOuterNormal(); | |
106 | 648456 | preGradX *= dim == dimWorld ? 0.5*a : 0.5*sqrt(a); | |
107 | preGradX /= preGradX.two_norm2(); | ||
108 | |||
109 |
2/2✓ Branch 0 taken 1296912 times.
✓ Branch 1 taken 648456 times.
|
1945368 | for (int compIdx = 0; compIdx < numComponents; compIdx++) |
110 | { | ||
111 | if constexpr (!FluidSystem::isTracerFluidSystem()) | ||
112 |
2/2✓ Branch 0 taken 648456 times.
✓ Branch 1 taken 648456 times.
|
1296912 | if (compIdx == FluidSystem::getMainComponent(phaseIdx)) |
113 | 648456 | continue; | |
114 | |||
115 | // interpolate mole fraction to scvf integration point | ||
116 | 648456 | Scalar x = 0.0; | |
117 |
4/4✓ Branch 0 taken 2593824 times.
✓ Branch 1 taken 648456 times.
✓ Branch 2 taken 2593824 times.
✓ Branch 3 taken 648456 times.
|
6484560 | for (const auto& scv : scvs(fvGeometry)) |
118 | 5187648 | x += massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx)*shapeValues[scv.indexInElement()][0]; | |
119 | |||
120 | // compute the diffusive flux by means of a finite difference | ||
121 | 648456 | auto gradX = preGradX; | |
122 | 648456 | gradX *= (massOrMoleFraction(facetVolVars, referenceSystem, phaseIdx, compIdx) - x); | |
123 | |||
124 | 648456 | componentFlux[compIdx] = -1.0*rho*Extrusion::area(fvGeometry, scvf) | |
125 | 648456 | *insideVolVars.extrusionFactor() | |
126 | 2593824 | *vtmv(scvf.unitOuterNormal(), | |
127 | facetVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx), | ||
128 | gradX); | ||
129 | |||
130 | if constexpr (!FluidSystem::isTracerFluidSystem()) | ||
131 | 648456 | if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx)) | |
132 | 1945368 | componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx]; | |
133 | } | ||
134 | |||
135 | 648456 | return componentFlux; | |
136 | } | ||
137 | |||
138 | // on interior Dirichlet boundaries use the facet mass/mole fraction and evaluate flux | ||
139 | ✗ | else if (bcTypes.hasOnlyDirichlet()) | |
140 | { | ||
141 | ✗ | for (int compIdx = 0; compIdx < numComponents; compIdx++) | |
142 | { | ||
143 | if constexpr (!FluidSystem::isTracerFluidSystem()) | ||
144 | ✗ | if (compIdx == FluidSystem::getMainComponent(phaseIdx)) | |
145 | ✗ | continue; | |
146 | |||
147 | // create vector with nodal mole/mass fractions | ||
148 | ✗ | std::vector<Scalar> xFractions(element.subEntities(dim)); | |
149 | ✗ | for (const auto& scv : scvs(fvGeometry)) | |
150 | ✗ | xFractions[scv.localDofIndex()] = massOrMoleFraction(elemVolVars[scv], referenceSystem, phaseIdx, compIdx); | |
151 | |||
152 | // substitute with facet mole/mass fractions for those scvs touching this facet | ||
153 | ✗ | for (const auto& scvfJ : scvfs(fvGeometry)) | |
154 | ✗ | if (scvfJ.interiorBoundary() && scvfJ.facetIndexInElement() == scvf.facetIndexInElement()) | |
155 | ✗ | xFractions[ fvGeometry.scv(scvfJ.insideScvIdx()).localDofIndex() ] | |
156 | ✗ | = massOrMoleFraction(problem.couplingManager().getLowDimVolVars(element, scvfJ), referenceSystem, phaseIdx, compIdx); | |
157 | |||
158 | // evaluate gradX at integration point | ||
159 | ✗ | Dune::FieldVector<Scalar, dimWorld> gradX(0.0); | |
160 | ✗ | for (const auto& scv : scvs(fvGeometry)) | |
161 | ✗ | gradX.axpy(xFractions[scv.localDofIndex()], fluxVarCache.gradN(scv.indexInElement())); | |
162 | |||
163 | // apply matrix diffusion coefficient and return the flux | ||
164 | ✗ | componentFlux[compIdx] = -1.0*rho*Extrusion::area(fvGeometry, scvf) | |
165 | ✗ | *insideVolVars.extrusionFactor() | |
166 | ✗ | *vtmv(scvf.unitOuterNormal(), | |
167 | insideVolVars.effectiveDiffusionCoefficient(phaseIdx, phaseIdx, compIdx), | ||
168 | gradX); | ||
169 | |||
170 | if constexpr (!FluidSystem::isTracerFluidSystem()) | ||
171 | ✗ | if (BalanceEqOpts::mainComponentIsBalanced(phaseIdx)) | |
172 | ✗ | componentFlux[FluidSystem::getMainComponent(phaseIdx)] -= componentFlux[compIdx]; | |
173 | } | ||
174 | |||
175 | ✗ | return componentFlux; | |
176 | } | ||
177 | |||
178 | // mixed boundary types are not supported | ||
179 | else | ||
180 | ✗ | DUNE_THROW(Dune::NotImplemented, "Mixed boundary types are not supported"); | |
181 | } | ||
182 | |||
183 | // compute transmissibilities ti for analytical jacobians | ||
184 | template<class Problem, class ElementVolumeVariables, class FluxVarCache> | ||
185 | static std::vector<Scalar> calculateTransmissibilities(const Problem& problem, | ||
186 | const Element& element, | ||
187 | const FVElementGeometry& fvGeometry, | ||
188 | const ElementVolumeVariables& elemVolVars, | ||
189 | const SubControlVolumeFace& scvf, | ||
190 | const FluxVarCache& fluxVarCache, | ||
191 | unsigned int phaseIdx) | ||
192 | { | ||
193 | DUNE_THROW(Dune::NotImplemented, "transmissibilty computation for BoxFacetCouplingFicksLaw"); | ||
194 | } | ||
195 | }; | ||
196 | |||
197 | } // end namespace Dumux | ||
198 | |||
199 | #endif | ||
200 |