3.2-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
discretization/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
30namespace Dumux {
31
36template<class BTypes>
38{
39public:
40 using BoundaryTypes = BTypes;
41
50 template<class Problem, class Element, class FVElementGeometry>
51 void update(const Problem &problem,
52 const Element &element,
53 const FVElementGeometry &fvGeometry)
54 {
55 using GridGeometry = typename FVElementGeometry::GridGeometry;
56 using GridView = typename GridGeometry::GridView;
57
58 vertexBCTypes_.resize( element.subEntities(GridView::dimension) );
59
60 hasDirichlet_ = false;
61 hasNeumann_ = false;
62 hasOutflow_ = false;
63
64 for (const auto& scv : scvs(fvGeometry))
65 {
66 int scvIdxLocal = scv.localDofIndex();
67 vertexBCTypes_[scvIdxLocal].reset();
68
69 if (fvGeometry.gridGeometry().dofOnBoundary(scv.dofIndex()))
70 {
71 vertexBCTypes_[scvIdxLocal] = problem.boundaryTypes(element, scv);
72
73 hasDirichlet_ = hasDirichlet_ || vertexBCTypes_[scvIdxLocal].hasDirichlet();
74 hasNeumann_ = hasNeumann_ || vertexBCTypes_[scvIdxLocal].hasNeumann();
75 hasOutflow_ = hasOutflow_ || vertexBCTypes_[scvIdxLocal].hasOutflow();
76 }
77 }
78 }
79
84 bool hasDirichlet() const
85 { return hasDirichlet_; }
86
91 bool hasNeumann() const
92 { return hasNeumann_; }
93
98 bool hasOutflow() const
99 { return hasOutflow_; }
100
101 /*
102 * \brief Access operator
103 * \return BoundaryTypes
104 */
105 const BoundaryTypes& operator[] (std::size_t i) const
106 {
107 assert(i < vertexBCTypes_.size());
108 return vertexBCTypes_[i];
109 }
110
111protected:
112 std::vector< BoundaryTypes > vertexBCTypes_;
113 bool hasDirichlet_ = false;
114 bool hasNeumann_ = false;
115 bool hasOutflow_ = false;
116};
117
118} // namespace Dumux
119
120#endif
Definition: adapt.hh:29
This class stores an array of BoundaryTypes objects.
Definition: discretization/box/elementboundarytypes.hh:38
std::vector< BoundaryTypes > vertexBCTypes_
Definition: discretization/box/elementboundarytypes.hh:112
bool hasOutflow() const
Returns whether the element potentially features an outflow boundary segment.
Definition: discretization/box/elementboundarytypes.hh:98
BTypes BoundaryTypes
Definition: discretization/box/elementboundarytypes.hh:40
bool hasDirichlet_
Definition: discretization/box/elementboundarytypes.hh:113
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry)
Update the boundary types for all vertices of an element.
Definition: discretization/box/elementboundarytypes.hh:51
const BoundaryTypes & operator[](std::size_t i) const
Definition: discretization/box/elementboundarytypes.hh:105
bool hasDirichlet() const
Returns whether the element has a vertex which contains a Dirichlet value.
Definition: discretization/box/elementboundarytypes.hh:84
bool hasNeumann_
Definition: discretization/box/elementboundarytypes.hh:114
bool hasNeumann() const
Returns whether the element potentially features a Neumann boundary segment.
Definition: discretization/box/elementboundarytypes.hh:91
bool hasOutflow_
Definition: discretization/box/elementboundarytypes.hh:115