3.4
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
box/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_BOX_ELEMENT_BOUNDARY_TYPES_HH
25#define DUMUX_BOX_ELEMENT_BOUNDARY_TYPES_HH
26
27#include <cassert>
28#include <vector>
29#include <dune/common/deprecated.hh>
30
31namespace Dumux {
32
37template<class BTypes>
39{
40public:
41 using BoundaryTypes = BTypes;
42
51 template<class Problem, class Element, class FVElementGeometry>
52 void update(const Problem &problem,
53 const Element &element,
54 const FVElementGeometry &fvGeometry)
55 {
56 using GridGeometry = typename FVElementGeometry::GridGeometry;
57 using GridView = typename GridGeometry::GridView;
58
59 vertexBCTypes_.resize( element.subEntities(GridView::dimension) );
60
61 hasDirichlet_ = false;
62 hasNeumann_ = false;
63 hasOutflow_ = false;
64
65 for (const auto& scv : scvs(fvGeometry))
66 {
67 int scvIdxLocal = scv.localDofIndex();
68 vertexBCTypes_[scvIdxLocal].reset();
69
70 if (fvGeometry.gridGeometry().dofOnBoundary(scv.dofIndex()))
71 {
72 vertexBCTypes_[scvIdxLocal] = problem.boundaryTypes(element, scv);
73
74 hasDirichlet_ = hasDirichlet_ || vertexBCTypes_[scvIdxLocal].hasDirichlet();
75 hasNeumann_ = hasNeumann_ || vertexBCTypes_[scvIdxLocal].hasNeumann();
76DUNE_NO_DEPRECATED_BEGIN
77 hasOutflow_ = hasOutflow_ || vertexBCTypes_[scvIdxLocal].hasOutflow();
78DUNE_NO_DEPRECATED_END
79 }
80 }
81 }
82
87 bool hasDirichlet() const
88 { return hasDirichlet_; }
89
94 bool hasNeumann() const
95 { return hasNeumann_; }
96
101 [[deprecated("Will be removed after release 3.4")]]
102 bool hasOutflow() const
103 { return hasOutflow_; }
104
105 /*
106 * \brief Access operator
107 * \return BoundaryTypes
108 */
109 const BoundaryTypes& operator[] (std::size_t i) const
110 {
111 assert(i < vertexBCTypes_.size());
112 return vertexBCTypes_[i];
113 }
114
115protected:
116 std::vector< BoundaryTypes > vertexBCTypes_;
117 bool hasDirichlet_ = false;
118 bool hasNeumann_ = false;
119 bool hasOutflow_ = false;
120};
121
122} // namespace Dumux
123
124#endif
Definition: adapt.hh:29
This class stores an array of BoundaryTypes objects.
Definition: box/elementboundarytypes.hh:39
std::vector< BoundaryTypes > vertexBCTypes_
Definition: box/elementboundarytypes.hh:116
bool hasOutflow() const
Returns whether the element potentially features an outflow boundary segment.
Definition: box/elementboundarytypes.hh:102
BTypes BoundaryTypes
Definition: box/elementboundarytypes.hh:41
bool hasDirichlet_
Definition: box/elementboundarytypes.hh:117
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry)
Update the boundary types for all vertices of an element.
Definition: box/elementboundarytypes.hh:52
const BoundaryTypes & operator[](std::size_t i) const
Definition: box/elementboundarytypes.hh:109
bool hasDirichlet() const
Returns whether the element has a vertex which contains a Dirichlet value.
Definition: box/elementboundarytypes.hh:87
bool hasNeumann_
Definition: box/elementboundarytypes.hh:118
bool hasNeumann() const
Returns whether the element potentially features a Neumann boundary segment.
Definition: box/elementboundarytypes.hh:94
bool hasOutflow_
Definition: box/elementboundarytypes.hh:119