24#ifndef DUMUX_PYTHON_DISCRETIZATION_GRIDGEOMETRY_HH
25#define DUMUX_PYTHON_DISCRETIZATION_GRIDGEOMETRY_HH
28#include <dune/common/classname.hh>
29#include <dune/python/pybind11/pybind11.h>
30#include <dune/python/common/typeregistry.hh>
37void registerSubControlVolume(pybind11::handle scope)
39 using namespace Dune::Python;
41 auto [cls, addedToRegistry] = insertClass<SCV>(
42 scope,
"SubControlVolume",
43 GenerateTypeName(Dune::className<SCV>()),
44 IncludeFiles{
"dumux/python/discretization/gridgeometry.hh"}
49 using pybind11::operator
""_a;
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);
62void registerSubControlVolumeFace(pybind11::handle scope)
64 using namespace Dune::Python;
66 auto [cls, addedToRegistry] = insertClass<SCVF>(
67 scope,
"SubControlVolumeFace",
68 GenerateTypeName(Dune::className<SCVF>()),
69 IncludeFiles{
"dumux/python/discretization/gridgeometry.hh"}
74 using pybind11::operator
""_a;
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);
88void registerFVElementGeometry(pybind11::handle scope)
90 using namespace Dune::Python;
92 auto [cls, addedToRegistry] = insertClass<FVEG>(
93 scope,
"FVElementGeometry",
94 GenerateTypeName(Dune::className<FVEG>()),
95 IncludeFiles{
"dumux/python/discretization/gridgeometry.hh"}
100 using pybind11::operator
""_a;
102 cls.def_property_readonly(
"numScvf", &FVEG::numScv);
103 cls.def_property_readonly(
"numScv", &FVEG::numScv);
104 cls.def_property_readonly(
"hasBoundaryScvf", &FVEG::hasBoundaryScvf);
106 using Element =
typename FVEG::Element;
107 cls.def(
"bind", [](FVEG& self,
const Element& element){
110 cls.def(
"bindElement", [](FVEG& self,
const Element& element){
111 self.bindElement(element);
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>());
128template <
class GG,
class... Options>
131 using pybind11::operator
""_a;
133 using GridView =
typename GG::GridView;
135 cls.def(pybind11::init([](
const GridView& gridView){
136 return std::make_shared<GG>(gridView);
139 cls.def(
"update", [](GG& self,
const GridView& gridView){
140 return self.update(gridView);
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);
150 cls.def_property_readonly_static(
"discMethod", [](
const pybind11::object&){
151 return GG::discMethod.name();
154 cls.def_property_readonly(
"localView", [](GG& self){
158 using LocalView =
typename GG::LocalView;
159 using Element =
typename LocalView::Element;
160 cls.def(
"boundLocalView", [](GG& self,
const Element& element){
166 using SubControlVolume =
typename GG::SubControlVolume;
167 Impl::registerSubControlVolume<SubControlVolume>(scope);
169 using SubControlVolumeFace =
typename GG::SubControlVolumeFace;
170 Impl::registerSubControlVolumeFace<SubControlVolumeFace>(scope);
173 Impl::registerFVElementGeometry<LocalView>(scope);
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