3.1-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
multidomain/facet/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 *****************************************************************************/
25#ifndef DUMUX_FACET_COUPLING_BOX_ELEMENT_BOUNDARY_TYPES_HH
26#define DUMUX_FACET_COUPLING_BOX_ELEMENT_BOUNDARY_TYPES_HH
27
28#include <cassert>
29#include <vector>
30
32
33namespace Dumux {
34
41template<class BTypes>
43{
44
45public:
57 template<class Problem, class Element, class FVElementGeometry>
58 void update(const Problem &problem,
59 const Element &element,
60 const FVElementGeometry &fvGeometry)
61 {
62 using GridGeometry = typename FVElementGeometry::GridGeometry;
63 using GridView = typename GridGeometry::GridView;
64
65 this->vertexBCTypes_.resize( element.subEntities(GridView::dimension) );
66
67 this->hasDirichlet_ = false;
68 this->hasNeumann_ = false;
69 this->hasOutflow_ = false;
70
71 for (const auto& scv : scvs(fvGeometry))
72 {
73 int scvIdxLocal = scv.localDofIndex();
74 this->vertexBCTypes_[scvIdxLocal].reset();
75
76 // lambda to update the element boundary info
77 auto updateElemBCInfo = [&] (const BTypes& bcTypes)
78 {
79 this->hasDirichlet_ = this->hasDirichlet_ || this->vertexBCTypes_[scvIdxLocal].hasDirichlet();
80 this->hasNeumann_ = this->hasNeumann_ || this->vertexBCTypes_[scvIdxLocal].hasNeumann();
81 this->hasOutflow_ = this->hasOutflow_ || this->vertexBCTypes_[scvIdxLocal].hasOutflow();
82 };
83
84 // We have to have Neumann-type BCs on nodes that touch interior boundaries
85 if (fvGeometry.gridGeometry().dofOnInteriorBoundary(scv.dofIndex()))
86 {
87 this->vertexBCTypes_[scvIdxLocal].setAllNeumann();
88 updateElemBCInfo(this->vertexBCTypes_[scvIdxLocal]);
89 }
90
91 // otherwise, let the problem decide
92 else if (fvGeometry.gridGeometry().dofOnBoundary(scv.dofIndex()))
93 {
94 this->vertexBCTypes_[scvIdxLocal] = problem.boundaryTypes(element, scv);
95 updateElemBCInfo(this->vertexBCTypes_[scvIdxLocal]);
96 }
97 }
98 }
99};
100
101} // namespace Dumux
102
103#endif
make the local view function available whenever we use the grid geometry
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 hasDirichlet_
Definition: discretization/box/elementboundarytypes.hh:113
bool hasNeumann_
Definition: discretization/box/elementboundarytypes.hh:114
bool hasOutflow_
Definition: discretization/box/elementboundarytypes.hh:115
This class stores an array of BoundaryTypes objects on an element for the box scheme with coupling oc...
Definition: multidomain/facet/box/elementboundarytypes.hh:43
void update(const Problem &problem, const Element &element, const FVElementGeometry &fvGeometry)
Update the boundary types for all vertices of an element.
Definition: multidomain/facet/box/elementboundarytypes.hh:58
Boundary types gathered on an element.