Loading [MathJax]/extensions/tex2jax.js

dune-mmesh (1.4)

utility.hh
1#ifndef DUNE_PYTHON_MMESH_UTILITY
2#define DUNE_PYTHON_MMESH_UTILITY
3
4#if HAVE_DUNE_FEM
5
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>
13
14namespace Dune
15{
16
17 namespace Fem
18 {
19
21 // Interface indicator //
23
24 template <class BulkGV>
25 struct InterfaceIndicator
26 : public BindableGridFunctionWithSpace<FemPy::GridPart<BulkGV>, Dune::FieldVector<double, 1>>
27 {
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;
34
35 InterfaceIndicator(const BulkGV &bulkGV)
36 : Base(FemPy::gridPart<BulkGV>(bulkGV), "interfaceindicator", 0),
37 onInterface_(false)
38 {}
39
40 void bind(const typename Base::EntityType &entity)
41 {
42 Base::bind(entity);
43 }
44
45 void bind(const typename BulkGV::Intersection &intersection, IntersectionSide side)
46 {
47 if( Base::gridPart().grid().isInterface( intersection ) )
48 onInterface_ = true;
49 else
50 onInterface_ = false;
51 }
52
53 template <class Point>
54 void evaluate(const Point &x, typename Base::RangeType &ret) const
55 {
56 if (onInterface_)
57 ret = RangeType(1);
58 else
59 ret = RangeType(0);
60 }
61
62 template <class Point>
63 void jacobian(const Point &x, typename Base::JacobianRangeType &ret) const
64 {
65 ret = typename Base::JacobianRangeType(0);
66 }
67
68 template <class Point>
69 void hessian(const Point &x, typename Base::HessianRangeType &ret) const
70 {
71 ret = typename Base::HessianRangeType(0);
72 }
73
74 private:
75 bool onInterface_;
76 };
77
78 template< class BulkGV >
79 inline static void registerInterfaceIndicator(pybind11::module module, pybind11::class_< InterfaceIndicator<BulkGV> > cls)
80 {
81 using pybind11::operator""_a;
82
83 cls.def( pybind11::init( [] ( const BulkGV &bulkGV ) {
84 return InterfaceIndicator<BulkGV> ( bulkGV );
85 } ), "bulkGV"_a, pybind11::keep_alive< 1, 2 >() );
86
87 cls.def_property_readonly( "scalar", [] ( InterfaceIndicator<BulkGV> &self ) { return true; } );
88
89 Dune::FemPy::registerGridFunction( module, cls );
90 }
91
92
93
95 // Normals //
97
98 template <class InterfaceGV>
99 struct Normals
100 : public BindableGridFunctionWithSpace<FemPy::GridPart<InterfaceGV>, Dune::FieldVector<double, InterfaceGV::dimensionworld>>
101 {
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;
108
109 Normals(const InterfaceGV &interfaceGV)
110 : Base(FemPy::gridPart<InterfaceGV>(interfaceGV), "normals", 0)
111 {}
112
113 void bind(const typename Base::EntityType &entity)
114 {
115 Base::bind(entity);
116 normal_ = Base::gridPart().grid().getMMesh().asIntersection(entity).centerUnitOuterNormal();
117 }
118
119 public:
120 template <class Point>
121 void evaluate(const Point &x, RangeType &ret) const
122 {
123 ret = normal_;
124 }
125
126 template <class Point>
127 void jacobian(const Point &x, typename Base::JacobianRangeType &ret) const
128 {
129 ret = typename Base::JacobianRangeType(0);
130 }
131
132 template <class Point>
133 void hessian(const Point &x, typename Base::HessianRangeType &ret) const
134 {
135 ret = typename Base::HessianRangeType(0);
136 }
137
138 private:
139 Dune::FieldVector<double, dimw> normal_;
140 };
141
142 template< class InterfaceGV >
143 inline static void registerNormals(pybind11::module module, pybind11::class_< Normals<InterfaceGV> > cls)
144 {
145 using pybind11::operator""_a;
146
147 cls.def( pybind11::init( [] ( const InterfaceGV &interfaceGV ) {
148 return Normals<InterfaceGV> ( interfaceGV );
149 } ), "interfaceGV"_a, pybind11::keep_alive< 1, 2 >() );
150
151 cls.def_property_readonly( "scalar", [] ( Normals<InterfaceGV> &self ) { return false; } );
152
153 Dune::FemPy::registerGridFunction( module, cls );
154 }
155
156
157
159 // DistanceFunction //
161
162 template <class GV>
163 struct DistanceFunction
164 : public BindableGridFunctionWithSpace<FemPy::GridPart<GV>, Dune::FieldVector<double, 1>>
165 {
166 using GridView = GV;
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;
172
173 DistanceFunction(const GridView &gridView)
174 : Base(FemPy::gridPart<GridView>(gridView), "distance", 0),
175 interpolation_(GeometryTypes::simplex(dim), std::vector<Dune::FieldVector<double, 1>>())
176 {}
177
178 void bind(const typename Base::EntityType &entity)
179 {
180 Base::bind(entity);
181
182 const auto& distance = this->gridPart().grid().indicator().distance();
183
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 )
187 {
188 const auto& vertex = entity.template subEntity<dim>( i );
189 distances[i] = distance(vertex);
190 }
191
192 interpolation_ = Dune::MultiLinearGeometry<double, dim, 1> (entity.type(), distances);
193 }
194
195 public:
196 template <class Point>
197 void evaluate(const Point &x, RangeType &ret) const
198 {
199 auto xLocal = Dune::Fem::coordinate(x);
200 ret = interpolation_.global(xLocal);
201 }
202
203 template <class Point>
204 void jacobian(const Point &x, typename Base::JacobianRangeType &ret) const
205 {
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];
210 }
211
212 template <class Point>
213 void hessian(const Point &x, typename Base::HessianRangeType &ret) const
214 {
215 ret = typename Base::HessianRangeType(0);
216 }
217
218 private:
219 Dune::MultiLinearGeometry<double, dim, 1> interpolation_;
220 };
221
222 template< class GridView >
223 inline static void registerDistanceFunction(pybind11::module module, pybind11::class_< DistanceFunction<GridView> > cls)
224 {
225 using pybind11::operator""_a;
226
227 cls.def( pybind11::init( [] ( const GridView &gridView ) {
228 return DistanceFunction<GridView> ( gridView );
229 } ), "gridView"_a, pybind11::keep_alive< 1, 2 >() );
230
231 cls.def_property_readonly( "scalar", [] ( DistanceFunction<GridView> &self ) { return true; } );
232
233 Dune::FemPy::registerGridFunction( module, cls );
234 }
235
236 } // namespace Fem
237
238} // namespace Dune
239
240#endif // HAVE_DUNE_FEM
241
242#endif // DUNE_PYTHON_MMESH_UTILITY
STL namespace.
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 5, 23:02, 2025)