3.6-git
DUNE for Multi-{Phase, Component, Scale, Physics, ...} flow and transport in porous media
python/discretization/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_property_readonly("center", &SCV::center);
52 cls.def_property_readonly("volume", &SCV::volume);
53 cls.def_property_readonly("dofIndex", &SCV::dofIndex);
54 cls.def_property_readonly("localDofIndex", &SCV::localDofIndex);
55 cls.def_property_readonly("dofPosition", &SCV::dofPosition);
56 cls.def_property_readonly("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_property_readonly("center", &SCVF::center);
77 cls.def_property_readonly("area", &SCVF::area);
78 cls.def_property_readonly("ipGlobal", &SCVF::ipGlobal);
79 cls.def_property_readonly("boundary", &SCVF::boundary);
80 cls.def_property_readonly("unitOuterNormal", &SCVF::unitOuterNormal);
81 cls.def_property_readonly("insideScvIdx", &SCVF::insideScvIdx);
82 cls.def_property_readonly("outsideScvIdx", [](SCVF& self){ return self.outsideScvIdx(); });
83 cls.def_property_readonly("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_property_readonly("numScvf", &FVEG::numScv);
103 cls.def_property_readonly("numScv", &FVEG::numScv);
104 cls.def_property_readonly("hasBoundaryScvf", &FVEG::hasBoundaryScvf);
105
106 using Element = typename FVEG::Element;
107 cls.def("bind", [](FVEG& self, const Element& element){
108 self.bind(element);
109 }, "element"_a);
110 cls.def("bindElement", [](FVEG& self, const Element& element){
111 self.bindElement(element);
112 }, "element"_a);
113
114 cls.def_property_readonly("scvs", [](FVEG& self){
115 const auto range = scvs(self);
116 return pybind11::make_iterator(range.begin(), range.end());
117 }, pybind11::keep_alive<0, 1>());
118 cls.def_property_readonly("scvfs", [](FVEG& self){
119 const auto range = scvfs(self);
120 return pybind11::make_iterator(range.begin(), range.end());
121 }, pybind11::keep_alive<0, 1>());
122 }
123}
124
125} // end namespace Impl
126
127// see python/dumux/discretization/__init__.py for how this is used for JIT compilation
128template <class GG, class... Options>
129void registerGridGeometry(pybind11::handle scope, pybind11::class_<GG, Options...> cls)
130{
131 using pybind11::operator""_a;
132
133 using GridView = typename GG::GridView;
134
135 cls.def(pybind11::init([](const GridView& gridView){
136 return std::make_shared<GG>(gridView);
137 }), "gridView"_a);
138
139 cls.def("update", [](GG& self, const GridView& gridView){
140 return self.update(gridView);
141 }, "gridView"_a);
142
143 cls.def_property_readonly("numDofs", &GG::numDofs);
144 cls.def_property_readonly("numScv", &GG::numScv);
145 cls.def_property_readonly("numScvf", &GG::numScvf);
146 cls.def_property_readonly("bBoxMax", &GG::bBoxMax);
147 cls.def_property_readonly("bBoxMin", &GG::bBoxMin);
148 cls.def_property_readonly("gridView", &GG::gridView);
149
150 cls.def_property_readonly_static("discMethod", [](const pybind11::object&){
151 return GG::discMethod.name();
152 });
153
154 cls.def_property_readonly("localView", [](GG& self){
155 return localView(self);
156 });
157
158 using LocalView = typename GG::LocalView;
159 using Element = typename LocalView::Element;
160 cls.def("boundLocalView", [](GG& self, const Element& element){
161 auto view = localView(self);
162 view.bind(element);
163 return view;
164 }, "element"_a);
165
166 using SubControlVolume = typename GG::SubControlVolume;
167 Impl::registerSubControlVolume<SubControlVolume>(scope);
168
169 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
170 Impl::registerSubControlVolumeFace<SubControlVolumeFace>(scope);
171
172 // also compile the corresponding local view
173 Impl::registerFVElementGeometry<LocalView>(scope);
174
175}
176
177} // namespace Dumux::Python
178
179#endif
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:171
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:36
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:38
Definition: python/assembly/fvassembler.hh:30
void registerGridGeometry(pybind11::handle scope, pybind11::class_< GG, Options... > cls)
Definition: python/discretization/gridgeometry.hh:129