Loading [MathJax]/extensions/tex2jax.js

dune-mmesh (1.4)

grid.hh
1// -*- tab-width: 2; indent-tabs-mode: nil; c-basic-offset: 2 -*-
2// vi: set et ts=4 sw=2 sts=2:
3#ifndef DUNE_PYTHON_MMESH_INTERFACEGRID_HH
4#define DUNE_PYTHON_MMESH_INTERFACEGRID_HH
5
6#include <array>
7#include <functional>
8#include <list>
9#include <map>
10#include <memory>
11#include <sstream>
12#include <type_traits>
13
14#include <dune/common/hybridutilities.hh>
15#include <dune/common/iteratorrange.hh>
16
17#include <dune/geometry/referenceelements.hh>
18#include <dune/geometry/type.hh>
19
20#include <dune/python/common/mpihelper.hh>
21#include <dune/python/common/typeregistry.hh>
22
23#include <dune/python/grid/capabilities.hh>
24#include <dune/python/grid/enums.hh>
25#include <dune/python/grid/factory.hh>
26#include <dune/python/grid/gridview.hh>
27#include <dune/python/grid/idset.hh>
28
29#include <dune/python/pybind11/functional.h>
30#include <dune/python/pybind11/numpy.h>
31#include <dune/python/pybind11/pybind11.h>
32#include <dune/python/pybind11/stl.h>
33
34namespace Dune
35{
36
37 namespace Python
38 {
39
40 namespace MMIFGrid
41 {
43 template< class Grid, class... options >
44 void registerHierarchicalGrid ( pybind11::module module, pybind11::class_<Grid, options... > cls )
45 {
46 Dune::Python::registerHierarchicalGrid( module, cls );
47
48 cls.def_property_readonly( "bulkGrid", [] ( const Grid &grid ) {
49 return grid.getMMesh().leafGridView();
50 },
51 R"doc(
52 Obtain bulk grid of the MMesh
53
54 Returns: bulk grid
55 )doc" );
56 }
57
58 } // end namespace MMIFGrid
59
60 namespace MMGrid
61 {
63 template< int d, class... options >
64 void registerHierarchicalGrid ( pybind11::module module, pybind11::class_<Dune::MovingMesh<d>, options...> cls )
65 {
66 Dune::Python::registerHierarchicalGrid( module, cls );
67
68 using Grid = Dune::MovingMesh< d >;
69
70 cls.def_property_readonly( "interfaceHierarchicalGrid", [] ( const Grid &grid ) -> const auto & {
71 return grid.interfaceGrid();
72 }, pybind11::return_value_policy::reference, pybind11::keep_alive< 0, 1 >(),
73 R"doc(
74 Obtain interface hierarchical grid of the MMesh
75
76 Returns: interface hierarchical grid
77 )doc" );
78
79 using Element = typename Grid::template Codim< 0 >::Entity;
80 using Vertex = typename Grid::template Codim< d >::Entity;
81 using Intersection = typename Grid::Intersection;
82 using InterfaceEntity = typename Grid::InterfaceEntity;
83 using InterfaceGrid = typename Grid::InterfaceGrid;
84 using InterfaceVertex = typename InterfaceGrid::template Codim< InterfaceGrid::dimension >::Entity;
85 using FieldVector = Dune::FieldVector< double, d >;
86
87 cls.def( "preAdapt", [] ( Grid &self ) {
88 self.preAdapt();
89 },
90 R"doc(
91 Prepare grid for adaption
92 )doc" );
93
94 cls.def( "ensureInterfaceMovement", [] ( Grid &self, const std::vector< FieldVector >& shifts ) {
95 return self.ensureInterfaceMovement( shifts );
96 },
97 R"doc(
98 Ensure the non-degeneration of the mesh after movement of the interface vertices
99 )doc" );
100
101 cls.def( "markElements", [] ( Grid &self ) {
102 return self.markElements();
103 },
104 R"doc(
105 Mark all elements in accordance to the default indicator
106 )doc" );
107
108 cls.def( "adapt", [] ( Grid &self ) {
109 self.adapt();
110 },
111 R"doc(
112 Adapt the grid
113 )doc" );
114
115 cls.def( "getConnectedComponent", [] ( Grid &self, const Element& element ) {
116 return self.getConnectedComponent( element );
117 },
118 R"doc(
119 Return the connected component of an entity
120 )doc" );
121
122 cls.def( "moveInterface", [] ( Grid &self, const std::vector< FieldVector >& shifts ) {
123 self.moveInterface( shifts );
124 },
125 R"doc(
126 Move the interface by the given movement for each interface vertex
127 )doc" );
128
129 cls.def( "moveVertices", [] ( Grid &self, const std::vector< FieldVector >& shifts ) {
130 self.moveVertices( shifts );
131 },
132 R"doc(
133 Move all vertices of the triangulation by the given movement
134 )doc" );
135
136 cls.def( "isInterface", [] ( Grid &self, const Intersection& intersection ) {
137 return self.isInterface( intersection );
138 },
139 R"doc(
140 Return if intersection is part of the interface
141 )doc" );
142
143 cls.def( "isInterface", [] ( Grid &self, const Vertex& vertex ) {
144 return self.isInterface( vertex );
145 },
146 R"doc(
147 Return if vertex is part of the interface
148 )doc" );
149
150 cls.def( "asInterfaceEntity", [] ( Grid &self, const Intersection& intersection ) {
151 return self.asInterfaceEntity( intersection );
152 },
153 R"doc(
154 Return intersection as entity of the interface grid
155 )doc" );
156
157 cls.def( "asIntersection", [] ( Grid &self, const InterfaceEntity& interfaceEntity ) {
158 return self.asIntersection( interfaceEntity );
159 },
160 R"doc(
161 Return entity of the interface grid as (some) intersection of the MMesh
162 )doc" );
163
164 cls.def( "addInterface", [] ( Grid &self, const Intersection& intersection ) {
165 self.addInterface( intersection );
166 },
167 R"doc(
168 Add the intersection to the set of interface edges
169 )doc" );
170
171 cls.def( "addInterface", [] ( Grid &self, const Intersection& intersection, const std::size_t marker ) {
172 self.addInterface( intersection, marker );
173 },
174 R"doc(
175 Add the intersection to the set of interface edges and mark it with marker
176 )doc" );
177
178 cls.def( "postAdapt", [] ( Grid &self ) {
179 self.postAdapt();
180 },
181 R"doc(
182 Remove adaption markers and connected components
183 )doc" );
184
185 cls.def( "isTip", [] ( Grid &self, const InterfaceVertex& interfaceVertex ) {
186 return interfaceVertex.impl().isTip();
187 },
188 R"doc(
189 Return if interface vertex is a tip
190 )doc" );
191
192 cls.def( "refineEdge", [] ( Grid &self, const Element& element, const std::size_t edgeIndex, const double where ) {
193 return self.refineEdge(element, edgeIndex, where);
194 },
195 R"doc(
196 Refine edge manually
197 )doc" );
198
199 cls.def( "refineEdge", [] ( Grid &self, const Element& element, const std::size_t edgeIndex ) {
200 return self.refineEdge(element, edgeIndex, 0.5);
201 },
202 R"doc(
203 Refine edge manually
204 )doc" );
205
206 cls.def( "removeVertex", [] ( Grid &self, const Vertex& vertex ) {
207 return self.removeVertex(vertex);
208 },
209 R"doc(
210 Remove vertex manually
211 )doc" );
212
213 cls.def( "removeVertex", [] ( Grid &self, const InterfaceVertex& vertex ) {
214 return self.removeVertex(vertex);
215 },
216 R"doc(
217 Remove interface vertex manually
218 )doc" );
219
220 cls.def( "insertVertexInCell", [] ( Grid &self, const FieldVector& position ) {
221 return self.insertVertexInCell(position);
222 },
223 R"doc(
224 Insert vertex in cell manually
225 )doc" );
226
227 }
228
229 } // namespace MMGrid
230
231 } // namespace Python
232
233} // namespace Dune
234
235#endif
The MMesh class templatized by the CGAL host grid type and the dimension.
Definition: mmesh.hh:140
Creative Commons License   |  Legal Statements / Impressum  |  Hosted by TU Dresden & Uni Heidelberg  |  generated with Hugo v0.111.3 (Apr 13, 22:42, 2025)