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 TwoPModel | ||
10 | * \copydoc Dumux::TwoPScvSaturationReconstruction | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_2P_SCV_SATURATION_RECONSTRUCTION_HH | ||
14 | #define DUMUX_2P_SCV_SATURATION_RECONSTRUCTION_HH | ||
15 | |||
16 | #include <dumux/discretization/method.hh> | ||
17 | |||
18 | namespace Dumux { | ||
19 | |||
20 | /*! | ||
21 | * \ingroup TwoPModel | ||
22 | * \brief Class that computes the nonwetting saturation in an scv from the saturation | ||
23 | * at the global degree of freedom. | ||
24 | * | ||
25 | * This is only necessary in conjunction with the box scheme where the degrees of | ||
26 | * freedom lie on material interfaces. There the nonwetting phase saturation is | ||
27 | * generally discontinuous. | ||
28 | */ | ||
29 | template<class DiscretizationMethod, bool enableReconstruction> | ||
30 | class TwoPScvSaturationReconstruction | ||
31 | { | ||
32 | public: | ||
33 | /*! | ||
34 | * \brief Compute the nonwetting phase saturation in an scv | ||
35 | * | ||
36 | * \note In the default case, we don't reconstruct anything. We do | ||
37 | * Reconstruction is only done when using the box method | ||
38 | * and enableReconstruction = true. | ||
39 | * | ||
40 | * \param spatialParams Class encapsulating the spatial parameters | ||
41 | * \param element The finite element the scv is embedded in | ||
42 | * \param scv The sub-control volume for which the saturation is computed | ||
43 | * \param elemSol The solution at all dofs inside this element | ||
44 | * \param sn The nonwetting phase saturation at the global dof | ||
45 | */ | ||
46 | template<class SpatialParams, class Element, class Scv, class ElemSol> | ||
47 | static typename ElemSol::PrimaryVariables::value_type | ||
48 | ✗ | reconstructSn(const SpatialParams& spatialParams, | |
49 | const Element& element, | ||
50 | const Scv& scv, | ||
51 | const ElemSol& elemSol, | ||
52 | typename ElemSol::PrimaryVariables::value_type sn) | ||
53 | ✗ | { return sn; } | |
54 | }; | ||
55 | |||
56 | //! Specialization for the box scheme with the interface solver enabled | ||
57 | template<> | ||
58 | class TwoPScvSaturationReconstruction<DiscretizationMethods::Box, /*enableReconstruction*/true> | ||
59 | { | ||
60 | public: | ||
61 | /*! | ||
62 | * \brief Compute the nonwetting phase saturation in an scv | ||
63 | * | ||
64 | * \param spatialParams Class encapsulating the spatial parameters | ||
65 | * \param element The finite element the scv is embedded in | ||
66 | * \param scv The sub-control volume for which the saturation is computed | ||
67 | * \param elemSol The solution at all dofs inside this element | ||
68 | * \param sn The nonwetting phase saturation at the global dof | ||
69 | */ | ||
70 | template<class SpatialParams, class Element, class Scv, class ElemSol> | ||
71 | static typename ElemSol::PrimaryVariables::value_type | ||
72 | 37748104 | reconstructSn(const SpatialParams& spatialParams, | |
73 | const Element& element, | ||
74 | const Scv& scv, | ||
75 | const ElemSol& elemSol, | ||
76 | typename ElemSol::PrimaryVariables::value_type sn) | ||
77 | { | ||
78 | // if this dof doesn't lie on a material interface, simply return Sn | ||
79 |
2/2✓ Branch 0 taken 5756912 times.
✓ Branch 1 taken 31991192 times.
|
37748104 | const auto& materialInterfaces = spatialParams.materialInterfaces(); |
80 |
4/4✓ Branch 0 taken 5756912 times.
✓ Branch 1 taken 31991192 times.
✓ Branch 2 taken 5756912 times.
✓ Branch 3 taken 31991192 times.
|
75496208 | if (!materialInterfaces.isOnMaterialInterface(scv)) |
81 | return sn; | ||
82 | |||
83 | // compute capillary pressure using material parameters associated with the dof | ||
84 | 5756912 | const auto& interfacePcSw = materialInterfaces.pcSwAtDof(scv); | |
85 | 5756912 | const auto pc = interfacePcSw.pc(/*ww=*/1.0 - sn); | |
86 | |||
87 | // reconstruct by inverting the pc-sw curve | ||
88 |
2/2✓ Branch 0 taken 2167952 times.
✓ Branch 1 taken 3413386 times.
|
5756912 | const auto& pcSw = spatialParams.fluidMatrixInteraction(element, scv, elemSol).pcSwCurve(); |
89 |
2/2✓ Branch 0 taken 3392238 times.
✓ Branch 1 taken 2189100 times.
|
5756912 | const auto pcMin = pcSw.endPointPc(); |
90 | |||
91 |
3/4✓ Branch 0 taken 3392238 times.
✓ Branch 1 taken 2189100 times.
✗ Branch 2 not taken.
✓ Branch 3 taken 3392238 times.
|
5581338 | if (pc < pcMin && pcMin > 0.0) |
92 | return 0.0; | ||
93 | else | ||
94 | 2364674 | return 1.0 - pcSw.sw(pc); | |
95 | } | ||
96 | }; | ||
97 | |||
98 | } // end namespace Dumux | ||
99 | |||
100 | #endif | ||
101 |