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-FileCopyrightText: 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 Assembly | ||
10 | * \brief Volume variables deflection helper. | ||
11 | */ | ||
12 | |||
13 | #ifndef DUMUX_ASSEMBLY_VOLVARS_DEFLECTION_HELPER_HH | ||
14 | #define DUMUX_ASSEMBLY_VOLVARS_DEFLECTION_HELPER_HH | ||
15 | |||
16 | #include <type_traits> | ||
17 | #include <dune/common/reservedvector.hh> | ||
18 | |||
19 | #ifndef DOXYGEN | ||
20 | |||
21 | namespace Dumux::Detail { | ||
22 | |||
23 | template<class VolVarAccessor, class FVElementGeometry> | ||
24 | class VolVarsDeflectionHelper | ||
25 | { | ||
26 | static constexpr int maxNumscvs = FVElementGeometry::maxNumElementScvs; | ||
27 | |||
28 | using SubControlVolume = typename FVElementGeometry::GridGeometry::SubControlVolume; | ||
29 | using VolumeVariablesRef = std::invoke_result_t<VolVarAccessor, const SubControlVolume&>; | ||
30 | using VolumeVariables = std::decay_t<VolumeVariablesRef>; | ||
31 | static_assert(std::is_lvalue_reference_v<VolumeVariablesRef> | ||
32 | && !std::is_const_v<std::remove_reference_t<VolumeVariablesRef>>); | ||
33 | |||
34 | public: | ||
35 | 14205695 | VolVarsDeflectionHelper(VolVarAccessor&& accessor, | |
36 | const FVElementGeometry& fvGeometry, | ||
37 | bool deflectAllVolVars) | ||
38 | 14205695 | : deflectAll_(deflectAllVolVars) | |
39 | 14205695 | , accessor_(std::move(accessor)) | |
40 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 697104 times.
|
14205695 | , fvGeometry_(fvGeometry) |
41 | { | ||
42 |
2/2✓ Branch 0 taken 604200 times.
✓ Branch 1 taken 11645126 times.
|
14205695 | if (deflectAll_) |
43 |
2/2✓ Branch 0 taken 1208400 times.
✓ Branch 1 taken 604200 times.
|
1812600 | for (const auto& curScv : scvs(fvGeometry)) |
44 | 1208400 | origVolVars_.push_back(accessor_(curScv)); | |
45 | 14205695 | } | |
46 | |||
47 | 41898290 | void setCurrent(const SubControlVolume& scv) | |
48 | { | ||
49 |
4/6✓ Branch 0 taken 37717642 times.
✓ Branch 1 taken 1208400 times.
✓ Branch 2 taken 2886504 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 30832 times.
✗ Branch 5 not taken.
|
41898290 | if (!deflectAll_) |
50 | { | ||
51 | 40689890 | origVolVars_.clear(); | |
52 | 40689890 | origVolVars_.push_back(accessor_(scv)); | |
53 | } | ||
54 | 15272594 | } | |
55 | |||
56 | template<class ElementSolution, class Problem> | ||
57 | 101619560 | void deflect(const ElementSolution& elemSol, | |
58 | const SubControlVolume& scv, | ||
59 | const Problem& problem) | ||
60 | { | ||
61 |
2/2✓ Branch 0 taken 3625200 times.
✓ Branch 1 taken 90272900 times.
|
101619560 | if (deflectAll_) |
62 |
2/2✓ Branch 0 taken 7250400 times.
✓ Branch 1 taken 3625200 times.
|
10875600 | for (const auto& curScv : scvs(fvGeometry_)) |
63 | 7250400 | accessor_(curScv).update(elemSol, problem, fvGeometry_.element(), curScv); | |
64 | else | ||
65 | 97994360 | accessor_(scv).update(elemSol, problem, fvGeometry_.element(), scv); | |
66 | 101619560 | } | |
67 | |||
68 | 100273640 | void restore(const SubControlVolume& scv) | |
69 | { | ||
70 |
2/2✓ Branch 0 taken 88926980 times.
✓ Branch 1 taken 3625200 times.
|
100273640 | if (!deflectAll_) |
71 | 96648440 | accessor_(scv) = origVolVars_[0]; | |
72 | else | ||
73 |
2/2✓ Branch 0 taken 7250400 times.
✓ Branch 1 taken 3625200 times.
|
10875600 | for (const auto& curScv : scvs(fvGeometry_)) |
74 | 7250400 | accessor_(curScv) = origVolVars_[curScv.localDofIndex()]; | |
75 | 100273640 | } | |
76 | |||
77 | private: | ||
78 | const bool deflectAll_; | ||
79 | VolVarAccessor accessor_; | ||
80 | const FVElementGeometry& fvGeometry_; | ||
81 | Dune::ReservedVector<VolumeVariables, maxNumscvs> origVolVars_; | ||
82 | }; | ||
83 | |||
84 | template<class Accessor, class FVElementGeometry> | ||
85 | VolVarsDeflectionHelper(Accessor&&, FVElementGeometry&&) -> VolVarsDeflectionHelper<Accessor, std::decay_t<FVElementGeometry>>; | ||
86 | |||
87 | } // end namespace Dumux::Detail | ||
88 | |||
89 | #endif // DOXYGEN | ||
90 | |||
91 | #endif | ||
92 |