3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
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 * See the file COPYING for full copying permissions. *
5 * *
6 * This program is free software: you can redistribute it and/or modify *
7 * it under the terms of the GNU General Public License as published by *
8 * the Free Software Foundation, either version 3 of the License, or *
9 * (at your option) any later version. *
10 * *
11 * This program is distributed in the hope that it will be useful, *
12 * but WITHOUT ANY WARRANTY; without even the implied warranty of *
13 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the *
14 * GNU General Public License for more details. *
15 * *
16 * You should have received a copy of the GNU General Public License *
17 * along with this program. If not, see <http://www.gnu.org/licenses/>. *
18 *****************************************************************************/
24#ifndef DUMUX_CVFE_ELEMENT_BOUNDARY_TYPES_HH
25#define DUMUX_CVFE_ELEMENT_BOUNDARY_TYPES_HH
26
27#include <cassert>
28#include <vector>
29
31
32namespace Dumux {
33
38template<class BTypes>
40{
41public:
42 using BoundaryTypes = BTypes;
43
52 template<class Problem, class FVElementGeometry>
53 void update(const Problem& problem,
54 const typename FVElementGeometry::Element& element,
55 const FVElementGeometry& fvGeometry)
56 {
57 using GridGeometry = typename FVElementGeometry::GridGeometry;
58
59 if constexpr (GridGeometry::discMethod == DiscretizationMethods::fcdiamond)
60 if (!fvGeometry.hasBoundaryScvf())
61 return;
62
63 bcTypes_.resize(fvGeometry.numScv());
64
65 hasDirichlet_ = false;
66 hasNeumann_ = false;
67
68 // boundary dofs always have boundary corresponding faces
69 if constexpr (GridGeometry::discMethod == DiscretizationMethods::fcdiamond)
70 {
71 for (const auto& scvf : scvfs(fvGeometry))
72 {
73 if (scvf.boundary())
74 {
75 const auto localIndex = fvGeometry.scv(scvf.insideScvIdx()).localDofIndex();
76 bcTypes_[localIndex] = problem.boundaryTypes(element, scvf);
77 hasDirichlet_ = hasDirichlet_ || bcTypes_[localIndex].hasDirichlet();
78 hasNeumann_ = hasNeumann_ || bcTypes_[localIndex].hasNeumann();
79 }
80 }
81 }
82
83 // generally, single vertices maybe on the boundary even if the face is not
84 // in this case we need a different mechanism
85 else
86 {
87 for (const auto& scv : scvs(fvGeometry))
88 {
89 const auto scvIdxLocal = scv.localDofIndex();
90 bcTypes_[scvIdxLocal].reset();
91
92 if (fvGeometry.gridGeometry().dofOnBoundary(scv.dofIndex()))
93 {
94 bcTypes_[scvIdxLocal] = problem.boundaryTypes(element, scv);
95 hasDirichlet_ = hasDirichlet_ || bcTypes_[scvIdxLocal].hasDirichlet();
96 hasNeumann_ = hasNeumann_ || bcTypes_[scvIdxLocal].hasNeumann();
97 }
98 }
99 }
100 }
101
106 bool hasDirichlet() const
107 { return hasDirichlet_; }
108
113 bool hasNeumann() const
114 { return hasNeumann_; }
115
116 /*
117 * \brief Access operator
118 * \return BoundaryTypes
119 */
120 [[deprecated("This operator is deprecated use get function")]]
121 const BoundaryTypes& operator[] (std::size_t i) const
122 {
123 assert(i < bcTypes_.size());
124 return bcTypes_[i];
125 }
126
127 /*
128 * \brief Access operator
129 * \return BoundaryTypes
130 * \note yields undefined behaviour of the scv is not on the boundary
131 */
132 template<class FVElementGeometry>
133 const BoundaryTypes& get(const FVElementGeometry& fvGeometry, const typename FVElementGeometry::SubControlVolumeFace& scvf) const
134 {
135 assert(scvf.boundary());
136 const auto localDofIdx = fvGeometry.scv(scvf.insideScvIdx()).localDofIndex();
137 assert(localDofIdx < bcTypes_.size());
138 return bcTypes_[localDofIdx];
139 }
140
141 /*
142 * \brief Access operator
143 * \return BoundaryTypes
144 * \note yields undefined behaviour of the scv is not on the boundary
145 */
146 template<class FVElementGeometry>
147 const BoundaryTypes& get(const FVElementGeometry&, const typename FVElementGeometry::SubControlVolume& scv) const
148 {
149 const auto localDofIdx = scv.localDofIndex();
150 assert(localDofIdx < bcTypes_.size());
151 return bcTypes_[localDofIdx];
152 }
153
154private:
155 std::vector<BoundaryTypes> bcTypes_;
156 bool hasDirichlet_ = false;
157 bool hasNeumann_ = false;
158};
159
160} // namespace Dumux
161
162#endif
The available discretization methods in Dumux.
Adaption of the non-isothermal two-phase two-component flow model to problems with CO2.
Definition: adapt.hh:29
constexpr FCDiamond fcdiamond
Definition: method.hh:141
This class stores an array of BoundaryTypes objects.
Definition: cvfe/elementboundarytypes.hh:40
BTypes BoundaryTypes
Definition: cvfe/elementboundarytypes.hh:42
const BoundaryTypes & operator[](std::size_t i) const
Definition: cvfe/elementboundarytypes.hh:121
bool hasDirichlet() const
Returns whether the element has a vertex which contains a Dirichlet value.
Definition: cvfe/elementboundarytypes.hh:106
const BoundaryTypes & get(const FVElementGeometry &, const typename FVElementGeometry::SubControlVolume &scv) const
Definition: cvfe/elementboundarytypes.hh:147
bool hasNeumann() const
Returns whether the element potentially features a Neumann boundary segment.
Definition: cvfe/elementboundarytypes.hh:113
const BoundaryTypes & get(const FVElementGeometry &fvGeometry, const typename FVElementGeometry::SubControlVolumeFace &scvf) const
Definition: cvfe/elementboundarytypes.hh:133
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:53