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 Typetraits | ||
10 | * \brief Type traits to be used for detecting new interfaces related to local dofs | ||
11 | */ | ||
12 | #ifndef DUMUX_TYPETRAITS_LOCALDOFS__HH | ||
13 | #define DUMUX_TYPETRAITS_LOCALDOFS__HH | ||
14 | |||
15 | #include <dune/common/std/type_traits.hh> | ||
16 | #include <dune/common/rangeutilities.hh> | ||
17 | |||
18 | #include <dumux/common/indextraits.hh> | ||
19 | #include <dumux/discretization/cvfe/localdof.hh> | ||
20 | |||
21 | namespace Dumux::Detail::LocalDofs { | ||
22 | |||
23 | //! helper struct detecting if a fvElementGeometry object has a numLocalDofs() function | ||
24 | template<class Imp> | ||
25 | using NumLocalDofsDetector = decltype( | ||
26 | std::declval<Imp>().numLocalDofs() | ||
27 | ); | ||
28 | |||
29 | template<typename FVElementGeometry> | ||
30 |
2/4✓ Branch 0 taken 2522976 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1092264 times.
✗ Branch 3 not taken.
|
21339236 | constexpr int numLocalDofs(const FVElementGeometry& fvGeometry) |
31 | { | ||
32 | if constexpr (Dune::Std::is_detected<NumLocalDofsDetector, FVElementGeometry>::value) | ||
33 | 997572 | return fvGeometry.numLocalDofs(); | |
34 | else | ||
35 | 14171801 | return fvGeometry.numScv(); | |
36 | } | ||
37 | |||
38 | //! helper struct detecting if a fvElementGeometry object has maxNumElementDofs | ||
39 | template<class Imp> | ||
40 | using MaxNumElementDofs = decltype( Imp::maxNumElementDofs ); | ||
41 | |||
42 | template<typename FVElementGeometry> | ||
43 | constexpr int maxNumLocalDofs() | ||
44 | { | ||
45 | if constexpr (Dune::Std::is_detected<MaxNumElementDofs, FVElementGeometry>::value) | ||
46 | return FVElementGeometry::maxNumElementDofs; | ||
47 | else | ||
48 | return FVElementGeometry::maxNumElementScvs; | ||
49 | } | ||
50 | |||
51 | //! helper struct detecting if a class has a localDofIndex() function | ||
52 | template<class Imp> | ||
53 | using LocalDofIndexDetector = decltype( | ||
54 | std::declval<Imp>().localDofIndex() | ||
55 | ); | ||
56 | |||
57 | template<class ScvOrLocalDof> | ||
58 |
2/4✓ Branch 2 taken 37712 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 18528 times.
✗ Branch 7 not taken.
|
17272452 | inline auto index(const ScvOrLocalDof& scvOrLocalDof) |
59 | { | ||
60 | if constexpr (Dune::Std::is_detected<LocalDofIndexDetector, ScvOrLocalDof>::value) | ||
61 | return scvOrLocalDof.localDofIndex(); | ||
62 | else | ||
63 |
2/4✓ Branch 2 taken 37712 times.
✗ Branch 3 not taken.
✓ Branch 6 taken 18528 times.
✗ Branch 7 not taken.
|
17272452 | return scvOrLocalDof.index(); |
64 | } | ||
65 | |||
66 | //! helper struct detecting if a fvElementGeometry object defines its own local dof type | ||
67 | template<class FVG> | ||
68 | using SpecifiesLocalDof = typename FVG::LocalDof; | ||
69 | |||
70 | template<class FVG> | ||
71 | using LocalDof_t = Dune::Std::detected_or_t< | ||
72 | Dumux::CVFE::LocalDof<typename IndexTraits<typename FVG::GridGeometry::GridView>::LocalIndex, | ||
73 | typename IndexTraits<typename FVG::GridGeometry::GridView>::GridIndex >, | ||
74 | SpecifiesLocalDof, | ||
75 | FVG | ||
76 | >; | ||
77 | |||
78 | } // end namespace Dumux::Detail::LocalDofs | ||
79 | |||
80 | #endif | ||
81 |