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::BoxMaterialInterfaces | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_2P_BOX_MATERIAL_INTERFACES_HH | ||
14 | #define DUMUX_2P_BOX_MATERIAL_INTERFACES_HH | ||
15 | |||
16 | #include <dune/common/exceptions.hh> | ||
17 | |||
18 | #include <dumux/discretization/method.hh> | ||
19 | #include <dumux/discretization/cvfe/elementsolution.hh> | ||
20 | |||
21 | namespace Dumux { | ||
22 | |||
23 | /*! | ||
24 | * \ingroup TwoPModel | ||
25 | * \brief Class that determines the material with the lowest capillary | ||
26 | * pressure (under fully water-saturated conditions) around the nodes | ||
27 | * of a grid. | ||
28 | * | ||
29 | * These parameters are then associated with the global degree | ||
30 | * of freedom. On the basis of these parameters, the saturations in the | ||
31 | * remaining sub-control volumes connected to the vertex can be reconstructed. | ||
32 | */ | ||
33 | template<class GridGeometry, class PcKrSw> | ||
34 | class BoxMaterialInterfaces | ||
35 | { | ||
36 | using SubControlVolume = typename GridGeometry::SubControlVolume; | ||
37 | |||
38 | public: | ||
39 | template<class SpatialParams, class SolutionVector> | ||
40 | 7 | BoxMaterialInterfaces(const GridGeometry& gridGeometry, | |
41 | const SpatialParams& spatialParams, | ||
42 | const SolutionVector& x) | ||
43 | 14 | { | |
44 |
1/2✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
|
7 | update(gridGeometry, spatialParams, x); |
45 | 7 | } | |
46 | |||
47 | /*! | ||
48 | * \brief Updates the scv -> dofparameter map | ||
49 | * | ||
50 | * \param gridGeometry The finite volume grid geometry | ||
51 | * \param spatialParams Class encapsulating the spatial parameters | ||
52 | * \param x The current state of the solution vector | ||
53 | */ | ||
54 | template<class SpatialParams, class SolutionVector> | ||
55 | 7 | void update(const GridGeometry& gridGeometry, | |
56 | const SpatialParams& spatialParams, | ||
57 | const SolutionVector& x) | ||
58 | { | ||
59 | // make sure this is only called for geometries of the box method! | ||
60 | if (GridGeometry::discMethod != DiscretizationMethods::box) | ||
61 | DUNE_THROW(Dune::InvalidStateException, "Determination of the interface material parameters with " | ||
62 | "this class only makes sense when using the box method!"); | ||
63 | |||
64 | 8 | isOnMaterialInterface_.resize(gridGeometry.numDofs(), false); | |
65 | 8 | pcSwAtDof_.resize(gridGeometry.numDofs(), nullptr); | |
66 |
1/2✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
|
10 | auto fvGeometry = localView(gridGeometry); |
67 |
14/18✓ Branch 1 taken 7 times.
✗ Branch 2 not taken.
✓ Branch 4 taken 7 times.
✗ Branch 5 not taken.
✓ Branch 6 taken 5963 times.
✓ Branch 7 taken 6 times.
✓ Branch 8 taken 4426 times.
✓ Branch 9 taken 3 times.
✓ Branch 10 taken 1539 times.
✓ Branch 11 taken 4426 times.
✓ Branch 12 taken 4426 times.
✓ Branch 13 taken 3 times.
✓ Branch 14 taken 8852 times.
✓ Branch 15 taken 3 times.
✓ Branch 17 taken 4426 times.
✗ Branch 18 not taken.
✓ Branch 20 taken 4426 times.
✗ Branch 21 not taken.
|
11942 | for (const auto& element : elements(gridGeometry.gridView())) |
68 | { | ||
69 |
1/2✓ Branch 1 taken 10388 times.
✗ Branch 2 not taken.
|
10388 | const auto elemSol = elementSolution(element, x, gridGeometry); |
70 |
1/2✓ Branch 1 taken 10388 times.
✗ Branch 2 not taken.
|
10388 | fvGeometry.bind(element); |
71 | |||
72 |
4/4✓ Branch 0 taken 39888 times.
✓ Branch 1 taken 10388 times.
✓ Branch 2 taken 39888 times.
✓ Branch 3 taken 10388 times.
|
60664 | for (const auto& scv : scvs(fvGeometry)) |
73 | { | ||
74 |
2/2✓ Branch 0 taken 2192 times.
✓ Branch 1 taken 31552 times.
|
39888 | const auto fluidMatrixInteraction = spatialParams.fluidMatrixInteraction(element, scv, elemSol); |
75 | 39888 | const auto& pcKrSw = fluidMatrixInteraction.pcSwCurve(); | |
76 | |||
77 | // assert current preconditions on requiring the spatial params to store the pckrSw curve | ||
78 | static_assert(std::is_lvalue_reference<typename std::decay_t<decltype(fluidMatrixInteraction)>::PcKrSwType>::value, | ||
79 | "In order to use the box-interface solver please provide access " | ||
80 | "to the material law parameters via returning (const) references"); | ||
81 | |||
82 | // if no parameters had been set, set them now | ||
83 |
4/4✓ Branch 0 taken 5497 times.
✓ Branch 1 taken 34391 times.
✓ Branch 2 taken 5497 times.
✓ Branch 3 taken 34391 times.
|
79776 | if (!pcSwAtDof_[scv.dofIndex()]) |
84 | 5497 | pcSwAtDof_[scv.dofIndex()] = &pcKrSw; | |
85 | |||
86 | // otherwise only use the current ones if endPointPc (e.g. Brooks-Corey entry pressure) is lower | ||
87 |
6/6✓ Branch 0 taken 486 times.
✓ Branch 1 taken 33905 times.
✓ Branch 2 taken 486 times.
✓ Branch 3 taken 33905 times.
✓ Branch 4 taken 360 times.
✓ Branch 5 taken 29504 times.
|
98646 | else if (pcKrSw.endPointPc() < pcSwAtDof_[scv.dofIndex()]->endPointPc()) |
88 | { | ||
89 | 360 | pcSwAtDof_[scv.dofIndex()] = &pcKrSw; | |
90 | 1080 | isOnMaterialInterface_[scv.dofIndex()] = true; | |
91 | } | ||
92 | |||
93 | // keep track of material interfaces in any case | ||
94 |
2/2✓ Branch 0 taken 3290 times.
✓ Branch 1 taken 30741 times.
|
34031 | else if ( !(pcKrSw == *(pcSwAtDof_[scv.dofIndex()])) ) |
95 | 9870 | isOnMaterialInterface_[scv.dofIndex()] = true; | |
96 | } | ||
97 | } | ||
98 | 7 | } | |
99 | |||
100 | //! Returns if this scv is connected to a material interface | ||
101 | bool isOnMaterialInterface(const SubControlVolume& scv) const | ||
102 |
4/4✓ Branch 0 taken 5756912 times.
✓ Branch 1 taken 31991192 times.
✓ Branch 2 taken 5756912 times.
✓ Branch 3 taken 31991192 times.
|
75496208 | { return isOnMaterialInterface_[scv.dofIndex()]; } |
103 | |||
104 | //! Returns the material parameters associated with the dof | ||
105 | const PcKrSw& pcSwAtDof(const SubControlVolume& scv) const | ||
106 | 11513824 | { return *(pcSwAtDof_[scv.dofIndex()]); } | |
107 | |||
108 | private: | ||
109 | std::vector<bool> isOnMaterialInterface_; | ||
110 | std::vector<const PcKrSw*> pcSwAtDof_; | ||
111 | }; | ||
112 | |||
113 | } // end namespace Dumux | ||
114 | |||
115 | #endif | ||
116 |