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 FaceCenteredStaggeredDiscretization | ||
10 | * \brief Boundary types gathered on an element | ||
11 | */ | ||
12 | #ifndef DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_BOUNDARY_TYPES_HH | ||
13 | #define DUMUX_DISCRETIZATION_FACECENTERED_STAGGERED_ELEMENT_BOUNDARY_TYPES_HH | ||
14 | |||
15 | #include <vector> | ||
16 | |||
17 | namespace Dumux { | ||
18 | |||
19 | /*! | ||
20 | * \ingroup FaceCenteredStaggeredDiscretization | ||
21 | * \brief This class stores an array of BoundaryTypes objects | ||
22 | */ | ||
23 | template<class BTypes> | ||
24 | 2761797 | class FaceCenteredStaggeredElementBoundaryTypes | |
25 | { | ||
26 | public: | ||
27 | using BoundaryTypes = BTypes; | ||
28 | |||
29 | /*! | ||
30 | * \brief Update the boundary types for all vertices of an element. | ||
31 | * | ||
32 | * \param problem The problem object which needs to be simulated | ||
33 | * \param element The DUNE Codim<0> entity for which the boundary | ||
34 | * types should be collected | ||
35 | * \param fvGeometry The element's finite volume geometry | ||
36 | */ | ||
37 | template<class Problem, class FVElementGeometry> | ||
38 | 2761797 | void update(const Problem& problem, | |
39 | const typename FVElementGeometry::Element& element, | ||
40 | const FVElementGeometry& fvGeometry) | ||
41 | { | ||
42 |
2/2✓ Branch 1 taken 171366 times.
✓ Branch 2 taken 2590431 times.
|
2761797 | if (!fvGeometry.hasBoundaryScvf()) |
43 | return; | ||
44 | |||
45 | 171366 | bcTypes_.resize(fvGeometry.numScvf()); | |
46 | |||
47 | 171366 | hasDirichlet_ = false; | |
48 | 171366 | hasNeumann_ = false; | |
49 | |||
50 |
4/4✓ Branch 0 taken 561154 times.
✓ Branch 1 taken 1782716 times.
✓ Branch 2 taken 2343870 times.
✓ Branch 3 taken 171366 times.
|
3064954 | for (const auto& scvf : scvfs(fvGeometry)) |
51 | { | ||
52 |
2/2✓ Branch 0 taken 561154 times.
✓ Branch 1 taken 1782716 times.
|
2343870 | if (scvf.boundary()) |
53 | { | ||
54 |
3/3✓ Branch 0 taken 160930 times.
✓ Branch 1 taken 366143 times.
✓ Branch 2 taken 34081 times.
|
561154 | bcTypes_[scvf.localIndex()] = problem.boundaryTypes(element, scvf); |
55 |
4/4✓ Branch 0 taken 221745 times.
✓ Branch 1 taken 339409 times.
✓ Branch 2 taken 148649 times.
✓ Branch 3 taken 73096 times.
|
561154 | hasDirichlet_ = hasDirichlet_ || bcTypes_[scvf.localIndex()].hasDirichlet(); |
56 |
4/4✓ Branch 0 taken 489910 times.
✓ Branch 1 taken 71244 times.
✓ Branch 2 taken 27790 times.
✓ Branch 3 taken 462120 times.
|
660188 | hasNeumann_ = hasNeumann_ || bcTypes_[scvf.localIndex()].hasNeumann(); |
57 | } | ||
58 | } | ||
59 | } | ||
60 | |||
61 | /*! | ||
62 | * \brief Returns whether the element has a vertex which contains | ||
63 | * a Dirichlet value. | ||
64 | */ | ||
65 | 319057751 | bool hasDirichlet() const | |
66 |
8/10✓ Branch 0 taken 355081 times.
✓ Branch 1 taken 711613 times.
✓ Branch 2 taken 18814298 times.
✓ Branch 3 taken 229557246 times.
✓ Branch 4 taken 2406022 times.
✓ Branch 5 taken 66048537 times.
✓ Branch 6 taken 40981 times.
✓ Branch 7 taken 1123973 times.
✗ Branch 8 not taken.
✗ Branch 9 not taken.
|
319057751 | { return hasDirichlet_; } |
67 | |||
68 | /*! | ||
69 | * \brief Returns whether the element potentially features a | ||
70 | * Neumann boundary segment. | ||
71 | */ | ||
72 | 315115526 | bool hasNeumann() const | |
73 |
5/6✗ Branch 0 not taken.
✓ Branch 1 taken 5589722 times.
✓ Branch 2 taken 902198 times.
✓ Branch 3 taken 13875044 times.
✓ Branch 4 taken 4687524 times.
✓ Branch 5 taken 290061038 times.
|
315115526 | { return hasNeumann_; } |
74 | |||
75 | /* | ||
76 | * \brief Access operator | ||
77 | * \return BoundaryTypes | ||
78 | */ | ||
79 | 10819837 | const BoundaryTypes& operator[] (std::size_t i) const | |
80 | { | ||
81 |
1/2✗ Branch 0 not taken.
✓ Branch 1 taken 10819837 times.
|
10819837 | assert(i < bcTypes_.size()); |
82 | 10819837 | return bcTypes_[i]; | |
83 | } | ||
84 | |||
85 | protected: | ||
86 | std::vector<BoundaryTypes> bcTypes_; | ||
87 | bool hasDirichlet_ = false; | ||
88 | bool hasNeumann_ = false; | ||
89 | }; | ||
90 | |||
91 | } // end namespace Dumux | ||
92 | |||
93 | #endif | ||
94 |