3.3.0
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
gridgeometry.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_PYTHON_DISCRETIZATION_GRIDGEOMETRY_HH
25#define DUMUX_PYTHON_DISCRETIZATION_GRIDGEOMETRY_HH
26
27#include <memory>
28#include <dune/common/classname.hh>
29#include <dune/python/pybind11/pybind11.h>
30#include <dune/python/common/typeregistry.hh>
31
32namespace Dumux::Python {
33
34namespace Impl {
35
36template<class SCV>
37void registerSubControlVolume(pybind11::handle scope)
38{
39 using namespace Dune::Python;
40
41 auto [cls, addedToRegistry] = insertClass<SCV>(
42 scope, "SubControlVolume",
43 GenerateTypeName(Dune::className<SCV>()),
44 IncludeFiles{"dumux/python/discretization/gridgeometry.hh"}
45 );
46
47 if (addedToRegistry)
48 {
49 using pybind11::operator""_a;
50
51 cls.def("center", &SCV::center);
52 cls.def("volume", &SCV::volume);
53 cls.def("dofIndex", &SCV::dofIndex);
54 cls.def("localDofIndex", &SCV::localDofIndex);
55 cls.def("dofPosition", &SCV::dofPosition);
56 cls.def("elementIndex", &SCV::elementIndex);
57
58 }
59}
60
61template<class SCVF>
62void registerSubControlVolumeFace(pybind11::handle scope)
63{
64 using namespace Dune::Python;
65
66 auto [cls, addedToRegistry] = insertClass<SCVF>(
67 scope, "SubControlVolumeFace",
68 GenerateTypeName(Dune::className<SCVF>()),
69 IncludeFiles{"dumux/python/discretization/gridgeometry.hh"}
70 );
71
72 if (addedToRegistry)
73 {
74 using pybind11::operator""_a;
75
76 cls.def("center", &SCVF::center);
77 cls.def("area", &SCVF::area);
78 cls.def("ipGlobal", &SCVF::ipGlobal);
79 cls.def("boundary", &SCVF::boundary);
80 cls.def("unitOuterNormal", &SCVF::unitOuterNormal);
81 cls.def("insideScvIdx", &SCVF::insideScvIdx);
82 cls.def("outsideScvIdx", [](SCVF& self){ return self.outsideScvIdx(); });
83 cls.def("index", &SCVF::index);
84 }
85}
86
87template<class FVEG>
88void registerFVElementGeometry(pybind11::handle scope)
89{
90 using namespace Dune::Python;
91
92 auto [cls, addedToRegistry] = insertClass<FVEG>(
93 scope, "FVElementGeometry",
94 GenerateTypeName(Dune::className<FVEG>()),
95 IncludeFiles{"dumux/python/discretization/gridgeometry.hh"}
96 );
97
98 if (addedToRegistry)
99 {
100 using pybind11::operator""_a;
101
102 cls.def("numScvf", &FVEG::numScv);
103 cls.def("numScv", &FVEG::numScv);
104 cls.def("bind", &FVEG::bind, "element"_a);
105 cls.def("hasBoundaryScvf", &FVEG::hasBoundaryScvf);
106 cls.def("scvs", [](FVEG& self){
107 const auto range = scvs(self);
108 return pybind11::make_iterator(range.begin(), range.end());
109 }, pybind11::keep_alive<0, 1>());
110 cls.def("scvfs", [](FVEG& self){
111 const auto range = scvfs(self);
112 return pybind11::make_iterator(range.begin(), range.end());
113 }, pybind11::keep_alive<0, 1>());
114 }
115}
116
117} // end namespace Impl
118
119// see python/dumux/discretization/__init__.py for how this is used for JIT compilation
120template <class GG, class... Options>
121void registerGridGeometry(pybind11::handle scope, pybind11::class_<GG, Options...> cls)
122{
123 using pybind11::operator""_a;
124
125 using GridView = typename GG::GridView;
126
127 cls.def(pybind11::init([](GridView gridView){
128 return std::make_shared<GG>(gridView);
129 }), "gridView"_a);
130
131 cls.def("update", &GG::update);
132 cls.def("numDofs", &GG::numDofs);
133 cls.def("numScv", &GG::numScv);
134 cls.def("numScvf", &GG::numScvf);
135
136 using SubControlVolume = typename GG::SubControlVolume;
137 Impl::registerSubControlVolume<SubControlVolume>(scope);
138
139 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
140 Impl::registerSubControlVolumeFace<SubControlVolumeFace>(scope);
141
142 // also compile the corresponding local view
143 using LocalView = typename GG::LocalView;
144 Impl::registerFVElementGeometry<LocalView>(scope);
145 // and make it accessible
146 cls.def("localView", [](GG& self){ return localView(self); });
147}
148
149} // namespace Dumux::Python
150
151#endif
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: python/common/boundarytypes.hh:33
void registerGridGeometry(pybind11::handle scope, pybind11::class_< GG, Options... > cls)
Definition: gridgeometry.hh:121