version 3.8
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// SPDX-FileCopyrightInfo: Copyright © DuMux Project contributors, see AUTHORS.md in root folder
5// SPDX-License-Identifier: GPL-3.0-or-later
6//
12#ifndef DUMUX_PYTHON_DISCRETIZATION_GRIDGEOMETRY_HH
13#define DUMUX_PYTHON_DISCRETIZATION_GRIDGEOMETRY_HH
14
15#include <memory>
16#include <dune/common/classname.hh>
17#include <dune/python/pybind11/pybind11.h>
18#include <dune/python/common/typeregistry.hh>
19
20namespace Dumux::Python {
21
22namespace Impl {
23
24template<class SCV>
25void registerSubControlVolume(pybind11::handle scope)
26{
27 using namespace Dune::Python;
28
29 auto [cls, addedToRegistry] = insertClass<SCV>(
30 scope, "SubControlVolume",
31 GenerateTypeName(Dune::className<SCV>()),
32 IncludeFiles{"dumux/python/discretization/gridgeometry.hh"}
33 );
34
35 if (addedToRegistry)
36 {
37 using pybind11::operator""_a;
38
39 cls.def_property_readonly("center", &SCV::center);
40 cls.def_property_readonly("volume", &SCV::volume);
41 cls.def_property_readonly("dofIndex", &SCV::dofIndex);
42 cls.def_property_readonly("localDofIndex", &SCV::localDofIndex);
43 cls.def_property_readonly("dofPosition", &SCV::dofPosition);
44 cls.def_property_readonly("elementIndex", &SCV::elementIndex);
45
46 }
47}
48
49template<class SCVF>
50void registerSubControlVolumeFace(pybind11::handle scope)
51{
52 using namespace Dune::Python;
53
54 auto [cls, addedToRegistry] = insertClass<SCVF>(
55 scope, "SubControlVolumeFace",
56 GenerateTypeName(Dune::className<SCVF>()),
57 IncludeFiles{"dumux/python/discretization/gridgeometry.hh"}
58 );
59
60 if (addedToRegistry)
61 {
62 using pybind11::operator""_a;
63
64 cls.def_property_readonly("center", &SCVF::center);
65 cls.def_property_readonly("area", &SCVF::area);
66 cls.def_property_readonly("ipGlobal", &SCVF::ipGlobal);
67 cls.def_property_readonly("boundary", &SCVF::boundary);
68 cls.def_property_readonly("unitOuterNormal", &SCVF::unitOuterNormal);
69 cls.def_property_readonly("insideScvIdx", &SCVF::insideScvIdx);
70 cls.def_property_readonly("outsideScvIdx", [](SCVF& self){ return self.outsideScvIdx(); });
71 cls.def_property_readonly("index", &SCVF::index);
72 }
73}
74
75template<class FVEG>
76void registerFVElementGeometry(pybind11::handle scope)
77{
78 using namespace Dune::Python;
79
80 auto [cls, addedToRegistry] = insertClass<FVEG>(
81 scope, "FVElementGeometry",
82 GenerateTypeName(Dune::className<FVEG>()),
83 IncludeFiles{"dumux/python/discretization/gridgeometry.hh"}
84 );
85
86 if (addedToRegistry)
87 {
88 using pybind11::operator""_a;
89
90 cls.def_property_readonly("numScvf", &FVEG::numScv);
91 cls.def_property_readonly("numScv", &FVEG::numScv);
92 cls.def_property_readonly("hasBoundaryScvf", &FVEG::hasBoundaryScvf);
93
94 using Element = typename FVEG::Element;
95 cls.def("bind", [](FVEG& self, const Element& element){
96 self.bind(element);
97 }, "element"_a);
98 cls.def("bindElement", [](FVEG& self, const Element& element){
99 self.bindElement(element);
100 }, "element"_a);
101
102 cls.def_property_readonly("scvs", [](FVEG& self){
103 const auto range = scvs(self);
104 return pybind11::make_iterator(range.begin(), range.end());
105 }, pybind11::keep_alive<0, 1>());
106 cls.def_property_readonly("scvfs", [](FVEG& self){
107 const auto range = scvfs(self);
108 return pybind11::make_iterator(range.begin(), range.end());
109 }, pybind11::keep_alive<0, 1>());
110 }
111}
112
113} // end namespace Impl
114
115// see python/dumux/discretization/__init__.py for how this is used for JIT compilation
116template <class GG, class... Options>
117void registerGridGeometry(pybind11::handle scope, pybind11::class_<GG, Options...> cls)
118{
119 using pybind11::operator""_a;
120
121 using GridView = typename GG::GridView;
122
123 cls.def(pybind11::init([](const GridView& gridView){
124 return std::make_shared<GG>(gridView);
125 }), "gridView"_a);
126
127 cls.def("update", [](GG& self, const GridView& gridView){
128 return self.update(gridView);
129 }, "gridView"_a);
130
131 cls.def_property_readonly("numDofs", &GG::numDofs);
132 cls.def_property_readonly("numScv", &GG::numScv);
133 cls.def_property_readonly("numScvf", &GG::numScvf);
134 cls.def_property_readonly("bBoxMax", &GG::bBoxMax);
135 cls.def_property_readonly("bBoxMin", &GG::bBoxMin);
136 cls.def_property_readonly("gridView", &GG::gridView);
137
138 cls.def_property_readonly_static("discMethod", [](const pybind11::object&){
139 return GG::discMethod.name();
140 });
141
142 cls.def_property_readonly("localView", [](GG& self){
143 return localView(self);
144 });
145
146 using LocalView = typename GG::LocalView;
147 using Element = typename LocalView::Element;
148 cls.def("boundLocalView", [](GG& self, const Element& element){
149 auto view = localView(self);
150 view.bind(element);
151 return view;
152 }, "element"_a);
153
154 using SubControlVolume = typename GG::SubControlVolume;
155 Impl::registerSubControlVolume<SubControlVolume>(scope);
156
157 using SubControlVolumeFace = typename GG::SubControlVolumeFace;
158 Impl::registerSubControlVolumeFace<SubControlVolumeFace>(scope);
159
160 // also compile the corresponding local view
161 Impl::registerFVElementGeometry<LocalView>(scope);
162
163}
164
165} // namespace Dumux::Python
166
167#endif
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition: localview.hh:26
auto volume(const Geometry &geo, unsigned int integrationOrder=4)
The volume of a given geometry.
Definition: volume.hh:159
Corners::value_type center(const Corners &corners)
The center of a given list of corners.
Definition: center.hh:24
Definition: python/assembly/fvassembler.hh:18
void registerGridGeometry(pybind11::handle scope, pybind11::class_< GG, Options... > cls)
Definition: python/discretization/gridgeometry.hh:117