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 CCMpfaDiscretization | ||
10 | * \brief Classes for sub control entities of the mpfa-o method. | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH | ||
13 | #define DUMUX_DISCRETIZATION_CC_MPFA_O_LOCAL_SUBCONTROLENTITIES_HH | ||
14 | |||
15 | #include <dune/common/fvector.hh> | ||
16 | |||
17 | namespace Dumux | ||
18 | { | ||
19 | |||
20 | /*! | ||
21 | * \ingroup CCMpfaDiscretization | ||
22 | * \brief Class for the interaction volume-local sub-control volume used | ||
23 | * in the mpfa-o scheme. | ||
24 | * | ||
25 | * \tparam IvIndexSet The type used for index sets within interaction volumes | ||
26 | * \tparam dim The dimensionality of the grid | ||
27 | * \tparam dimWorld The dimension of the world the grid is embedded in | ||
28 | */ | ||
29 | template< class IvIndexSet, class Scalar, int dim, int dimWorld> | ||
30 | class CCMpfaOInteractionVolumeLocalScv | ||
31 | { | ||
32 | |||
33 | public: | ||
34 | // export some types | ||
35 | using GridIndexType = typename IvIndexSet::GridIndexType; | ||
36 | using LocalIndexType = typename IvIndexSet::LocalIndexType; | ||
37 | using GlobalCoordinate = Dune::FieldVector<Scalar, dimWorld>; | ||
38 | using ctype = typename GlobalCoordinate::value_type; | ||
39 | using LocalBasis = std::array< GlobalCoordinate, dim >; | ||
40 | |||
41 | static constexpr int myDimension = dim; | ||
42 | static constexpr int worldDimension = dimWorld; | ||
43 | |||
44 | //! The default constructor | ||
45 | CCMpfaOInteractionVolumeLocalScv() = default; | ||
46 | |||
47 | /*! | ||
48 | * \brief The constructor | ||
49 | * | ||
50 | * \param helper Helper class for mpfa schemes | ||
51 | * \param fvGeometry The element finite volume geometry | ||
52 | * \param scv The grid sub-control volume | ||
53 | * \param localIndex The iv-local index of this scvIdx | ||
54 | * \param indexSet The interaction volume index set | ||
55 | */ | ||
56 | template<class MpfaHelper, class FVElementGeometry, class SubControlVolume> | ||
57 | 37370591 | CCMpfaOInteractionVolumeLocalScv(const MpfaHelper& helper, | |
58 | const FVElementGeometry& fvGeometry, | ||
59 | const SubControlVolume& scv, | ||
60 | const LocalIndexType localIndex, | ||
61 | const IvIndexSet& indexSet) | ||
62 | : indexSet_(&indexSet) | ||
63 | , globalScvIndex_(scv.dofIndex()) | ||
64 | 74741182 | , localDofIndex_(localIndex) | |
65 | { | ||
66 | // center of the global scv | ||
67 | 37370591 | const auto& center = scv.center(); | |
68 | |||
69 | // set up local basis | ||
70 | 37370591 | LocalBasis localBasis; | |
71 |
2/2✓ Branch 0 taken 74741182 times.
✓ Branch 1 taken 37370591 times.
|
112111773 | for (unsigned int coordIdx = 0; coordIdx < myDimension; ++coordIdx) |
72 | { | ||
73 | 74741182 | const auto scvfIdx = indexSet.nodalIndexSet().gridScvfIndex(localDofIndex_, coordIdx); | |
74 | 74741182 | const auto& scvf = fvGeometry.scvf(scvfIdx); | |
75 | 74741182 | localBasis[coordIdx] = scvf.ipGlobal(); | |
76 | 224223546 | localBasis[coordIdx] -= center; | |
77 | } | ||
78 | |||
79 | 37370591 | nus_ = helper.calculateInnerNormals(localBasis); | |
80 | 72033912 | detX_ = helper.calculateDetX(localBasis); | |
81 | 37370591 | } | |
82 | |||
83 | //! detX is needed for setting up the omegas in the interaction volumes | ||
84 | ✗ | ctype detX() const | |
85 | ✗ | { return detX_; } | |
86 | |||
87 | //! grid index related to this scv | ||
88 | ✗ | GridIndexType gridScvIndex() const | |
89 | ✗ | { return globalScvIndex_; } | |
90 | |||
91 | //! returns the index in the set of cell unknowns of the iv | ||
92 | ✗ | LocalIndexType localDofIndex() const | |
93 | ✗ | { return localDofIndex_; } | |
94 | |||
95 | //! iv-local index of the coordir's scvf in this scv | ||
96 | ✗ | LocalIndexType localScvfIndex(unsigned int coordDir) const | |
97 | { | ||
98 | ✗ | assert(coordDir < myDimension); | |
99 | 800 | return indexSet_->localScvfIndex(localDofIndex_, coordDir); | |
100 | } | ||
101 | |||
102 | //! the nu vectors are needed for setting up the omegas of the iv | ||
103 | const GlobalCoordinate& nu(unsigned int coordDir) const | ||
104 | { | ||
105 | assert(coordDir < myDimension); | ||
106 | 1631659208 | return nus_[coordDir]; | |
107 | } | ||
108 | |||
109 | private: | ||
110 | const IvIndexSet* indexSet_; | ||
111 | GridIndexType globalScvIndex_; | ||
112 | LocalIndexType localDofIndex_; | ||
113 | LocalBasis nus_; | ||
114 | ctype detX_; | ||
115 | }; | ||
116 | |||
117 | /*! | ||
118 | * \ingroup CCMpfaDiscretization | ||
119 | * \brief Class for the interaction volume-local sub-control volume face | ||
120 | * used in the mpfa-o scheme. | ||
121 | * | ||
122 | * \tparam IvIndexSet The type used for index sets within interaction volumes | ||
123 | */ | ||
124 | template< class IvIndexSet > | ||
125 | struct CCMpfaOInteractionVolumeLocalScvf | ||
126 | { | ||
127 | using ScvfNeighborLocalIndexSet = typename IvIndexSet::ScvfNeighborLocalIndexSet; | ||
128 | |||
129 | public: | ||
130 | // export index types | ||
131 | using GridIndexType = typename IvIndexSet::GridIndexType; | ||
132 | using LocalIndexType = typename IvIndexSet::LocalIndexType; | ||
133 | |||
134 | //! The default constructor | ||
135 | CCMpfaOInteractionVolumeLocalScvf() = default; | ||
136 | |||
137 | /*! | ||
138 | * \brief The constructor | ||
139 | * | ||
140 | * \param scvf The grid sub-control volume face | ||
141 | * \param localScvIndices The iv-local neighboring scv indices | ||
142 | * \param localDofIdx This scvf's interaction volume-local dof index | ||
143 | * \param isDirichlet Specifies if this scv is on a Dirichlet boundary | ||
144 | */ | ||
145 | template< class SubControlVolumeFace > | ||
146 | 39881298 | CCMpfaOInteractionVolumeLocalScvf(const SubControlVolumeFace& scvf, | |
147 | const ScvfNeighborLocalIndexSet& localScvIndices, | ||
148 | const LocalIndexType localDofIdx, | ||
149 | const bool isDirichlet) | ||
150 | : isDirichlet_(isDirichlet) | ||
151 | , scvfIdxGlobal_(scvf.index()) | ||
152 | , localDofIndex_(localDofIdx) | ||
153 |
0/17✗ Branch 0 not taken.
✗ Branch 1 not taken.
✗ Branch 2 not taken.
✗ Branch 4 not taken.
✗ Branch 5 not taken.
✗ Branch 7 not taken.
✗ Branch 8 not taken.
✗ Branch 10 not taken.
✗ Branch 11 not taken.
✗ Branch 14 not taken.
✗ Branch 15 not taken.
✗ Branch 18 not taken.
✗ Branch 19 not taken.
✗ Branch 22 not taken.
✗ Branch 23 not taken.
✗ Branch 26 not taken.
✗ Branch 27 not taken.
|
40070358 | , neighborScvIndicesLocal_(&localScvIndices) |
154 | {} | ||
155 | |||
156 | //! This is either the iv-local index of the intermediate unknown (interior/Neumann face) | ||
157 | //! or the index of the Dirichlet boundary within the vol vars (Dirichlet faces) | ||
158 | ✗ | LocalIndexType localDofIndex() const { return localDofIndex_; } | |
159 | |||
160 | //! returns the grid view-global index of this scvf | ||
161 | ✗ | GridIndexType gridScvfIndex() const { return scvfIdxGlobal_; } | |
162 | |||
163 | //! Returns the local indices of the scvs neighboring this scvf | ||
164 | ✗ | const ScvfNeighborLocalIndexSet& neighboringLocalScvIndices() const { return *neighborScvIndicesLocal_; } | |
165 | |||
166 | //! states if this is scvf is on a Dirichlet boundary | ||
167 | ✗ | bool isDirichlet() const { return isDirichlet_; } | |
168 | |||
169 | private: | ||
170 | bool isDirichlet_; | ||
171 | GridIndexType scvfIdxGlobal_; | ||
172 | LocalIndexType localDofIndex_; | ||
173 | const ScvfNeighborLocalIndexSet* neighborScvIndicesLocal_; | ||
174 | }; | ||
175 | |||
176 | } // end namespace Dumux | ||
177 | |||
178 | #endif | ||
179 |