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;
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);
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;
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);
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(
"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>());
120template <
class GG,
class... Options>
123 using pybind11::operator
""_a;
125 using GridView =
typename GG::GridView;
127 cls.def(pybind11::init([](GridView gridView){
128 return std::make_shared<GG>(gridView);
131 cls.def(
"update", &GG::update);
132 cls.def(
"numDofs", &GG::numDofs);
133 cls.def(
"numScv", &GG::numScv);
134 cls.def(
"numScvf", &GG::numScvf);
136 using SubControlVolume =
typename GG::SubControlVolume;
137 Impl::registerSubControlVolume<SubControlVolume>(scope);
139 using SubControlVolumeFace =
typename GG::SubControlVolumeFace;
140 Impl::registerSubControlVolumeFace<SubControlVolumeFace>(scope);
143 using LocalView =
typename GG::LocalView;
144 Impl::registerFVElementGeometry<LocalView>(scope);
146 cls.def(
"localView", [](GG& self){
return localView(self); });
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