1#ifndef DUNE_PYTHON_MMESH_UTILITY
2#define DUNE_PYTHON_MMESH_UTILITY
6#include <dune/common/exceptions.hh>
7#include <dune/python/common/typeregistry.hh>
8#include <dune/fem/common/intersectionside.hh>
9#include <dune/fem/function/localfunction/const.hh>
10#include <dune/fem/function/localfunction/bindable.hh>
11#include <dune/fempy/py/grid/gridpart.hh>
12#include <dune/fempy/py/grid/function.hh>
24 template <
class BulkGV>
25 struct InterfaceIndicator
26 :
public BindableGridFunctionWithSpace<FemPy::GridPart<BulkGV>, Dune::FieldVector<double, 1>>
28 using GridView = BulkGV;
29 using GridPartType = FemPy::GridPart<BulkGV>;
30 using Base = BindableGridFunctionWithSpace<GridPartType, Dune::FieldVector<double, 1>>;
31 using SideGeometry =
typename GridPartType::GridType::Intersection::LocalGeometry::Implementation;
32 using RangeType =
typename Base::RangeType;
33 static constexpr bool scalar =
true;
35 InterfaceIndicator(
const BulkGV &bulkGV)
36 : Base(FemPy::gridPart<BulkGV>(bulkGV),
"interfaceindicator", 0),
40 void bind(
const typename Base::EntityType &entity)
45 void bind(
const typename BulkGV::Intersection &intersection, IntersectionSide side)
47 if( Base::gridPart().grid().isInterface( intersection ) )
53 template <
class Po
int>
54 void evaluate(
const Point &x,
typename Base::RangeType &ret)
const
62 template <
class Po
int>
63 void jacobian(
const Point &x,
typename Base::JacobianRangeType &ret)
const
65 ret =
typename Base::JacobianRangeType(0);
68 template <
class Po
int>
69 void hessian(
const Point &x,
typename Base::HessianRangeType &ret)
const
71 ret =
typename Base::HessianRangeType(0);
78 template<
class BulkGV >
79 inline static void registerInterfaceIndicator(pybind11::module module, pybind11::class_< InterfaceIndicator<BulkGV> > cls)
81 using pybind11::operator
""_a;
83 cls.def( pybind11::init( [] (
const BulkGV &bulkGV ) {
84 return InterfaceIndicator<BulkGV> ( bulkGV );
85 } ),
"bulkGV"_a, pybind11::keep_alive< 1, 2 >() );
87 cls.def_property_readonly(
"scalar", [] ( InterfaceIndicator<BulkGV> &self ) {
return true; } );
89 Dune::FemPy::registerGridFunction( module, cls );
98 template <
class InterfaceGV>
100 :
public BindableGridFunctionWithSpace<FemPy::GridPart<InterfaceGV>, Dune::FieldVector<double, InterfaceGV::dimensionworld>>
102 static constexpr int dimw = InterfaceGV::dimensionworld;
103 using GridView = InterfaceGV;
104 using GridPartType = FemPy::GridPart<InterfaceGV>;
105 using Base = BindableGridFunctionWithSpace<GridPartType, Dune::FieldVector<double, dimw>>;
106 using RangeType =
typename Base::RangeType;
107 static constexpr bool scalar =
false;
109 Normals(
const InterfaceGV &interfaceGV)
110 : Base(FemPy::gridPart<InterfaceGV>(interfaceGV),
"normals", 0)
113 void bind(
const typename Base::EntityType &entity)
116 normal_ = Base::gridPart().grid().getMMesh().asIntersection(entity).centerUnitOuterNormal();
120 template <
class Po
int>
121 void evaluate(
const Point &x, RangeType &ret)
const
126 template <
class Po
int>
127 void jacobian(
const Point &x,
typename Base::JacobianRangeType &ret)
const
129 ret =
typename Base::JacobianRangeType(0);
132 template <
class Po
int>
133 void hessian(
const Point &x,
typename Base::HessianRangeType &ret)
const
135 ret =
typename Base::HessianRangeType(0);
139 Dune::FieldVector<double, dimw> normal_;
142 template<
class InterfaceGV >
143 inline static void registerNormals(pybind11::module module, pybind11::class_< Normals<InterfaceGV> > cls)
145 using pybind11::operator
""_a;
147 cls.def( pybind11::init( [] (
const InterfaceGV &interfaceGV ) {
148 return Normals<InterfaceGV> ( interfaceGV );
149 } ),
"interfaceGV"_a, pybind11::keep_alive< 1, 2 >() );
151 cls.def_property_readonly(
"scalar", [] ( Normals<InterfaceGV> &self ) {
return false; } );
153 Dune::FemPy::registerGridFunction( module, cls );
163 struct DistanceFunction
164 :
public BindableGridFunctionWithSpace<FemPy::GridPart<GV>, Dune::FieldVector<double, 1>>
167 static constexpr int dim = GridView::dimensionworld;
168 using GridPartType = FemPy::GridPart<GridView>;
169 using Base = BindableGridFunctionWithSpace<GridPartType, Dune::FieldVector<double, 1>>;
170 using RangeType =
typename Base::RangeType;
171 static constexpr bool scalar =
true;
173 DistanceFunction(
const GridView &gridView)
174 : Base(FemPy::gridPart<GridView>(gridView),
"distance", 0),
175 interpolation_(GeometryTypes::simplex(dim),
std::vector<Dune::FieldVector<double, 1>>())
178 void bind(
const typename Base::EntityType &entity)
182 const auto& distance = this->gridPart().grid().indicator().distance();
184 std::vector<Dune::FieldVector<double, 1>> distances;
185 distances.resize(dim+1);
186 for ( std::size_t i = 0; i < entity.subEntities( dim ); ++i )
188 const auto& vertex = entity.template subEntity<dim>( i );
189 distances[i] = distance(vertex);
192 interpolation_ = Dune::MultiLinearGeometry<double, dim, 1> (entity.type(), distances);
196 template <
class Po
int>
197 void evaluate(
const Point &x, RangeType &ret)
const
199 auto xLocal = Dune::Fem::coordinate(x);
200 ret = interpolation_.global(xLocal);
203 template <
class Po
int>
204 void jacobian(
const Point &x,
typename Base::JacobianRangeType &ret)
const
206 auto xLocal = Dune::Fem::coordinate(x);
207 auto jacT = interpolation_.jacobianTransposed(xLocal);
208 for (std::size_t i = 0; i < dim; ++i)
209 ret[0][i] = jacT[i][0];
212 template <
class Po
int>
213 void hessian(
const Point &x,
typename Base::HessianRangeType &ret)
const
215 ret =
typename Base::HessianRangeType(0);
219 Dune::MultiLinearGeometry<double, dim, 1> interpolation_;
222 template<
class Gr
idView >
223 inline static void registerDistanceFunction(pybind11::module module, pybind11::class_< DistanceFunction<GridView> > cls)
225 using pybind11::operator
""_a;
227 cls.def( pybind11::init( [] (
const GridView &gridView ) {
228 return DistanceFunction<GridView> ( gridView );
229 } ),
"gridView"_a, pybind11::keep_alive< 1, 2 >() );
231 cls.def_property_readonly(
"scalar", [] ( DistanceFunction<GridView> &self ) {
return true; } );
233 Dune::FemPy::registerGridFunction( module, cls );