12#ifndef DUMUX_PYTHON_DISCRETIZATION_GRIDGEOMETRY_HH 
   13#define DUMUX_PYTHON_DISCRETIZATION_GRIDGEOMETRY_HH 
   16#include <dune/common/classname.hh> 
   17#include <dune/python/pybind11/pybind11.h> 
   18#include <dune/python/common/typeregistry.hh> 
   25void registerSubControlVolume(pybind11::handle scope)
 
   27    using namespace Dune::Python;
 
   29    auto [cls, addedToRegistry] = insertClass<SCV>(
 
   30        scope, 
"SubControlVolume",
 
   31        GenerateTypeName(Dune::className<SCV>()),
 
   32        IncludeFiles{
"dumux/python/discretization/gridgeometry.hh"}
 
   37        using pybind11::operator
""_a;
 
   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);
 
   50void registerSubControlVolumeFace(pybind11::handle scope)
 
   52    using namespace Dune::Python;
 
   54    auto [cls, addedToRegistry] = insertClass<SCVF>(
 
   55        scope, 
"SubControlVolumeFace",
 
   56        GenerateTypeName(Dune::className<SCVF>()),
 
   57        IncludeFiles{
"dumux/python/discretization/gridgeometry.hh"}
 
   62        using pybind11::operator
""_a;
 
   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);
 
   76void registerFVElementGeometry(pybind11::handle scope)
 
   78    using namespace Dune::Python;
 
   80    auto [cls, addedToRegistry] = insertClass<FVEG>(
 
   81        scope, 
"FVElementGeometry",
 
   82        GenerateTypeName(Dune::className<FVEG>()),
 
   83        IncludeFiles{
"dumux/python/discretization/gridgeometry.hh"}
 
   88        using pybind11::operator
""_a;
 
   90        cls.def_property_readonly(
"numScvf", &FVEG::numScv);
 
   91        cls.def_property_readonly(
"numScv", &FVEG::numScv);
 
   92        cls.def_property_readonly(
"hasBoundaryScvf", &FVEG::hasBoundaryScvf);
 
   94        using Element = 
typename FVEG::Element;
 
   95        cls.def(
"bind", [](FVEG& self, 
const Element& element){
 
   98        cls.def(
"bindElement", [](FVEG& self, 
const Element& element){
 
   99            self.bindElement(element);
 
  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>());
 
  116template <
class GG, 
class... Options>
 
  119    using pybind11::operator
""_a;
 
  121    using GridView = 
typename GG::GridView;
 
  123    cls.def(pybind11::init([](
const GridView& gridView){
 
  124        return std::make_shared<GG>(gridView);
 
  127    cls.def(
"update", [](GG& self, 
const GridView& gridView){
 
  128        return self.update(gridView);
 
  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);
 
  138    cls.def_property_readonly_static(
"discMethod", [](
const pybind11::object&){
 
  139        return GG::discMethod.name();
 
  142    cls.def_property_readonly(
"localView", [](GG& self){
 
  146    using LocalView = 
typename GG::LocalView;
 
  147    using Element = 
typename LocalView::Element;
 
  148    cls.def(
"boundLocalView", [](GG& self, 
const Element& element){
 
  154    using SubControlVolume = 
typename GG::SubControlVolume;
 
  155    Impl::registerSubControlVolume<SubControlVolume>(scope);
 
  157    using SubControlVolumeFace = 
typename GG::SubControlVolumeFace;
 
  158    Impl::registerSubControlVolumeFace<SubControlVolumeFace>(scope);
 
  161    Impl::registerFVElementGeometry<LocalView>(scope);
 
 
GridCache::LocalView localView(const GridCache &gridCache)
Free function to get the local view of a grid cache object.
Definition localview.hh:26
Definition python/assembly/fvassembler.hh:18
void registerGridGeometry(pybind11::handle scope, pybind11::class_< GG, Options... > cls)
Definition python/discretization/gridgeometry.hh:117