version 3.8
cvfe/elementboundarytypes.hh
Go to the documentation of this file.
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//
12#ifndef DUMUX_CVFE_ELEMENT_BOUNDARY_TYPES_HH
13#define DUMUX_CVFE_ELEMENT_BOUNDARY_TYPES_HH
14
15#include <cassert>
16#include <vector>
17
19
20namespace Dumux {
21
26template<class BTypes>
28{
29public:
30 using BoundaryTypes = BTypes;
31
40 template<class Problem, class FVElementGeometry>
41 void update(const Problem& problem,
42 const typename FVElementGeometry::Element& element,
43 const FVElementGeometry& fvGeometry)
44 {
45 bcTypes_.resize(fvGeometry.numScv());
46
47 hasDirichlet_ = false;
48 hasNeumann_ = false;
49
50 for (const auto& scv : scvs(fvGeometry))
51 {
52 const auto scvIdxLocal = scv.localDofIndex();
53 bcTypes_[scvIdxLocal].reset();
54
55 if (fvGeometry.gridGeometry().dofOnBoundary(scv.dofIndex()))
56 {
57 bcTypes_[scvIdxLocal] = problem.boundaryTypes(element, scv);
58 hasDirichlet_ = hasDirichlet_ || bcTypes_[scvIdxLocal].hasDirichlet();
59 hasNeumann_ = hasNeumann_ || bcTypes_[scvIdxLocal].hasNeumann();
60 }
61 }
62 }
63
68 bool hasDirichlet() const
69 { return hasDirichlet_; }
70
75 bool hasNeumann() const
76 { return hasNeumann_; }
77
78 /*
79 * \brief Access operator
80 * \return BoundaryTypes
81 * \note yields undefined behaviour of the scv is not on the boundary
82 */
83 template<class FVElementGeometry>
84 const BoundaryTypes& get(const FVElementGeometry& fvGeometry, const typename FVElementGeometry::SubControlVolumeFace& scvf) const
85 {
86 assert(scvf.boundary());
87 const auto localDofIdx = fvGeometry.scv(scvf.insideScvIdx()).localDofIndex();
88 assert(localDofIdx < bcTypes_.size());
89 return bcTypes_[localDofIdx];
90 }
91
92 /*
93 * \brief Access operator
94 * \return BoundaryTypes
95 * \note yields undefined behaviour of the scv is not on the boundary
96 */
97 template<class FVElementGeometry>
98 const BoundaryTypes& get(const FVElementGeometry&, const typename FVElementGeometry::SubControlVolume& scv) const
99 {
100 const auto localDofIdx = scv.localDofIndex();
101 assert(localDofIdx < bcTypes_.size());
102 return bcTypes_[localDofIdx];
103 }
104
105private:
106 std::vector<BoundaryTypes> bcTypes_;
107 bool hasDirichlet_ = false;
108 bool hasNeumann_ = false;
109};
110
111} // namespace Dumux
112
113#endif
This class stores an array of BoundaryTypes objects.
Definition: cvfe/elementboundarytypes.hh:28
BTypes BoundaryTypes
Definition: cvfe/elementboundarytypes.hh:30
bool hasDirichlet() const
Returns whether the element has a vertex which contains a Dirichlet value.
Definition: cvfe/elementboundarytypes.hh:68
const BoundaryTypes & get(const FVElementGeometry &, const typename FVElementGeometry::SubControlVolume &scv) const
Definition: cvfe/elementboundarytypes.hh:98
bool hasNeumann() const
Returns whether the element potentially features a Neumann boundary segment.
Definition: cvfe/elementboundarytypes.hh:75
const BoundaryTypes & get(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf) const
Definition: cvfe/elementboundarytypes.hh:84
void update(const Problem &problem, const typename FVElementGeometry::Element &element, const FVElementGeometry &fvGeometry)
Update the boundary types for all vertices of an element.
Definition: cvfe/elementboundarytypes.hh:41
The available discretization methods in Dumux.
Definition: adapt.hh:17