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 CVFEDiscretization | ||
10 | * \brief The local element solution class for control-volume finite element methods | ||
11 | */ | ||
12 | #ifndef DUMUX_CVFE_ELEMENT_SOLUTION_HH | ||
13 | #define DUMUX_CVFE_ELEMENT_SOLUTION_HH | ||
14 | |||
15 | #include <type_traits> | ||
16 | #include <dune/common/reservedvector.hh> | ||
17 | |||
18 | #include <dumux/common/typetraits/localdofs_.hh> | ||
19 | #include <dumux/discretization/method.hh> | ||
20 | #include <dumux/discretization/localdoftraits.hh> | ||
21 | #include <dumux/discretization/cvfe/localdof.hh> | ||
22 | |||
23 | namespace Dumux { | ||
24 | |||
25 | /*! | ||
26 | * \ingroup CVFEDiscretization | ||
27 | * \brief The element solution vector | ||
28 | */ | ||
29 | template<class FVElementGeometry, class PV> | ||
30 | class CVFEElementSolution | ||
31 | { | ||
32 | using GridGeometry = typename FVElementGeometry::GridGeometry; | ||
33 | using GridView = typename GridGeometry::GridView; | ||
34 | using Element = typename GridView::template Codim<0>::Entity; | ||
35 | |||
36 | static constexpr int dim = GridView::dimension; | ||
37 | |||
38 | // Dofs are located at corners and in the element center | ||
39 | static constexpr int numCubeDofs = Detail::LocalDofTraits< | ||
40 | GridView, typename GridGeometry::DiscretizationMethod | ||
41 | >::numCubeElementDofs; | ||
42 | |||
43 | public: | ||
44 | //! export the primary variables type | ||
45 | using PrimaryVariables = PV; | ||
46 | |||
47 | //! Default constructor | ||
48 | 9 | CVFEElementSolution() = default; | |
49 | |||
50 | //! Constructor with element and solution and grid geometry | ||
51 | template<class SolutionVector> | ||
52 | 86771357 | CVFEElementSolution(const Element& element, const SolutionVector& sol, | |
53 | const GridGeometry& gridGeometry) | ||
54 | 127730484 | { | |
55 |
10/14✓ Branch 2 taken 1782467 times.
✓ Branch 3 taken 24492 times.
✓ Branch 7 taken 34908 times.
✗ Branch 8 not taken.
✓ Branch 10 taken 69700 times.
✗ Branch 11 not taken.
✓ Branch 5 taken 30156 times.
✓ Branch 6 taken 20 times.
✓ Branch 4 taken 12372 times.
✓ Branch 13 taken 12279 times.
✗ Branch 14 not taken.
✓ Branch 16 taken 34 times.
✗ Branch 17 not taken.
✓ Branch 1 taken 346732 times.
|
86771257 | update(element, sol, gridGeometry); |
56 | 18790643 | } | |
57 | |||
58 | //! Constructor with element and elemVolVars and fvGeometry | ||
59 | template<class ElementVolumeVariables> | ||
60 | 10542417 | CVFEElementSolution(const Element& element, const ElementVolumeVariables& elemVolVars, | |
61 | const FVElementGeometry& fvGeometry) | ||
62 | 10542417 | { | |
63 |
2/4✓ Branch 0 taken 2522976 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1092264 times.
✗ Branch 3 not taken.
|
10542417 | priVars_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry)); |
64 |
3/4✓ Branch 0 taken 11890472 times.
✓ Branch 1 taken 4384440 times.
✓ Branch 2 taken 1092264 times.
✗ Branch 3 not taken.
|
19986807 | for (const auto& scv : scvs(fvGeometry)) |
65 | 9444390 | priVars_[scv.localDofIndex()] = elemVolVars[scv].priVars(); | |
66 | 4124844 | } | |
67 | |||
68 | //! extract the element solution from the solution vector using a mapper | ||
69 | template<class SolutionVector> | ||
70 | 130055725 | void update(const Element& element, const SolutionVector& sol, | |
71 | const GridGeometry& gridGeometry) | ||
72 | { | ||
73 | // TODO: this implementation only works if there is only one dof per codim/entity | ||
74 | // As local-global mappings are provided by the grid geometry | ||
75 | // this logic should maybe become part of the grid geometry. | ||
76 | 130055725 | const auto& localCoeff = gridGeometry.feCache().get(element.type()).localCoefficients(); | |
77 | 130055725 | priVars_.resize(localCoeff.size()); | |
78 |
3/3✓ Branch 0 taken 255478800 times.
✓ Branch 1 taken 135081124 times.
✓ Branch 2 taken 16132227 times.
|
594501635 | for (int localDofIdx = 0; localDofIdx < localCoeff.size(); ++localDofIdx) |
79 | { | ||
80 | 464445910 | const auto& localKey = localCoeff.localKey(localDofIdx); | |
81 | 464445910 | priVars_[localDofIdx] = sol[ | |
82 | 464445910 | gridGeometry.dofMapper().subIndex( | |
83 | 438758042 | element, localKey.subEntity(), localKey.codim() | |
84 | ) | ||
85 | ]; | ||
86 | } | ||
87 | 130055725 | } | |
88 | |||
89 | //! extract the element solution from the solution vector using a local fv geometry | ||
90 | template<class SolutionVector> | ||
91 | void update(const Element& element, const SolutionVector& sol, | ||
92 | const FVElementGeometry& fvGeometry) | ||
93 | { | ||
94 | priVars_.resize(Detail::LocalDofs::numLocalDofs(fvGeometry)); | ||
95 | for (const auto& localDof : localDofs(fvGeometry)) | ||
96 | priVars_[localDof.index()] = sol[localDof.dofIndex()]; | ||
97 | } | ||
98 | |||
99 | //! return the size of the element solution | ||
100 | 13311096 | std::size_t size() const | |
101 |
2/2✓ Branch 0 taken 6655548 times.
✓ Branch 1 taken 6655548 times.
|
13311096 | { return priVars_.size(); } |
102 | |||
103 | //! bracket operator const access | ||
104 | template<typename IndexType> | ||
105 | 924597130 | const PrimaryVariables& operator [](IndexType i) const | |
106 |
5/9✓ Branch 3 taken 2595772 times.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 6 not taken.
✓ Branch 1 taken 11358321 times.
✓ Branch 2 taken 36120 times.
✓ Branch 7 taken 93600 times.
✗ Branch 8 not taken.
✓ Branch 0 taken 800797 times.
|
922515592 | { return priVars_[i]; } |
107 | |||
108 | //! bracket operator access | ||
109 | template<typename IndexType> | ||
110 | 286863612 | PrimaryVariables& operator [](IndexType i) | |
111 |
15/21✓ Branch 0 taken 226808 times.
✓ Branch 1 taken 4903336 times.
✓ Branch 4 taken 4298382 times.
✓ Branch 5 taken 874560 times.
✓ Branch 10 taken 2509828 times.
✗ Branch 11 not taken.
✓ Branch 13 taken 2042736 times.
✗ Branch 14 not taken.
✓ Branch 19 taken 17472 times.
✗ Branch 20 not taken.
✓ Branch 22 taken 2940304 times.
✗ Branch 23 not taken.
✓ Branch 3 taken 68939250 times.
✓ Branch 8 taken 824256 times.
✓ Branch 9 taken 82944 times.
✓ Branch 7 taken 22382 times.
✓ Branch 16 taken 1056 times.
✓ Branch 17 taken 100000 times.
✗ Branch 18 not taken.
✓ Branch 12 taken 82944 times.
✗ Branch 2 not taken.
|
278327224 | { return priVars_[i]; } |
112 | |||
113 | private: | ||
114 | Dune::ReservedVector<PrimaryVariables, numCubeDofs> priVars_; | ||
115 | }; | ||
116 | |||
117 | //! Make an element solution for control-volume finite element schemes | ||
118 | template<class Element, class SolutionVector, class GridGeometry> | ||
119 |
4/6✓ Branch 2 taken 1688867 times.
✓ Branch 3 taken 24472 times.
✓ Branch 5 taken 30156 times.
✗ Branch 6 not taken.
✓ Branch 4 taken 12279 times.
✗ Branch 1 not taken.
|
85103781 | auto elementSolution(const Element& element, const SolutionVector& sol, const GridGeometry& gg) |
120 | -> std::enable_if_t<DiscretizationMethods::isCVFE<typename GridGeometry::DiscretizationMethod>, | ||
121 | CVFEElementSolution<typename GridGeometry::LocalView, | ||
122 | std::decay_t<decltype(std::declval<SolutionVector>()[0])>> | ||
123 | > | ||
124 | { | ||
125 | using PrimaryVariables = std::decay_t<decltype(std::declval<SolutionVector>()[0])>; | ||
126 |
11/15✓ Branch 3 taken 1781501 times.
✓ Branch 4 taken 34715 times.
✓ Branch 6 taken 73528 times.
✗ Branch 7 not taken.
✓ Branch 2 taken 2393474 times.
✓ Branch 5 taken 2438154 times.
✓ Branch 1 taken 404370 times.
✓ Branch 9 taken 44820 times.
✗ Branch 10 not taken.
✓ Branch 8 taken 12869 times.
✓ Branch 11 taken 1794158 times.
✗ Branch 12 not taken.
✓ Branch 14 taken 34 times.
✗ Branch 15 not taken.
✓ Branch 0 taken 534520 times.
|
87725613 | return CVFEElementSolution<typename GridGeometry::LocalView, PrimaryVariables>(element, sol, gg); |
127 | } | ||
128 | |||
129 | //! Make an element solution for control-volume finite element schemes | ||
130 | template<class Element, class ElementVolumeVariables, class FVElementGeometry> | ||
131 |
2/4✓ Branch 0 taken 2522976 times.
✗ Branch 1 not taken.
✓ Branch 2 taken 1092264 times.
✗ Branch 3 not taken.
|
10525533 | auto elementSolution(const Element& element, const ElementVolumeVariables& elemVolVars, const FVElementGeometry& gg) |
132 | -> std::enable_if_t<DiscretizationMethods::isCVFE<typename FVElementGeometry::GridGeometry::DiscretizationMethod>, | ||
133 | CVFEElementSolution<FVElementGeometry, | ||
134 | typename ElementVolumeVariables::VolumeVariables::PrimaryVariables>> | ||
135 | { | ||
136 | using PrimaryVariables = typename ElementVolumeVariables::VolumeVariables::PrimaryVariables; | ||
137 |
5/6✓ Branch 0 taken 2528731 times.
✓ Branch 1 taken 279356 times.
✓ Branch 2 taken 1092264 times.
✗ Branch 3 not taken.
✓ Branch 4 taken 5045952 times.
✓ Branch 5 taken 2522976 times.
|
15577239 | return CVFEElementSolution<FVElementGeometry, PrimaryVariables>(element, elemVolVars, gg); |
138 | } | ||
139 | |||
140 | } // end namespace Dumux | ||
141 | |||
142 | #endif | ||
143 |