File: referenceelements.hh

package info (click to toggle)
dune-geometry 2.11.0-1
  • links: PTS, VCS
  • area: main
  • in suites: forky
  • size: 2,316 kB
  • sloc: cpp: 15,164; python: 262; makefile: 6
file content (269 lines) | stat: -rw-r--r-- 11,191 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
// SPDX-FileCopyrightText: Copyright © DUNE Project contributors, see file LICENSE.md in module root
// SPDX-License-Identifier: LicenseRef-GPL-2.0-only-with-DUNE-exception
#ifndef DUNE_PYTHON_GEOMETRY_REFERENCEELEMENTS_HH
#define DUNE_PYTHON_GEOMETRY_REFERENCEELEMENTS_HH

#include <array>
#include <functional>
#include <string>

#include <dune/common/visibility.hh>

#include <dune/geometry/referenceelements.hh>
#include <dune/geometry/type.hh>

#include <dune/python/common/fvecmatregistry.hh>

#include <dune/python/geometry/quadraturerules.hh>
#include <dune/python/pybind11/pybind11.h>

namespace Dune
{

  namespace Python
  {

    namespace detail
    {

      // referenceElementSize
      // --------------------

      template< class RefElement >
      inline static int referenceElementSize ( const RefElement &refElement, int c )
      {
        if( (c < 0) || (c > RefElement::dimension) )
          throw pybind11::value_error( "Invalid codimension: " + std::to_string( c ) + " (must be in [0, " + std::to_string( RefElement::dimension ) + "])." );
        return refElement.size( c );
      }

      template< class RefElement >
      inline static int referenceElementSize ( const RefElement &refElement, int i, int c, int cc )
      {
        const int size = detail::referenceElementSize( refElement, c );
        if( (i < 0) || (i >= size) )
          throw pybind11::value_error( "Invalid index: " + std::to_string( i ) + " (must be in [0, " + std::to_string( size ) + "))." );
        if( (cc < c) || (cc > RefElement::dimension) )
          throw pybind11::value_error( "Invalid codimension: " + std::to_string( cc ) + " (must be in [" + std::to_string( c ) + ", " + std::to_string( RefElement::dimension ) + "])." );
        return refElement.size( i, c, cc );
      }



      // referenceElementSubEntity
      // -------------------------

      template< class RefElement >
      inline static int referenceElementSubEntity ( const RefElement &refElement, int i, int c, int ii, int cc )
      {
        const int size = detail::referenceElementSize( refElement, i, c, cc );
        if( (ii < 0) || (ii >= size) )
          throw pybind11::value_error( "Invalid index: " + std::to_string( i ) + " (must be in [0, " + std::to_string( size ) + "))." );
        return refElement.subEntity( i, c, ii, cc );
      }



      // referenceElementPosition
      // ------------------------

      template< class RefElement >
      inline static auto referenceElementPosition ( const RefElement &refElement, int i, int c )
      {
        const int size = detail::referenceElementSize( refElement, c );
        if( (i < 0) || (i >= size) )
          throw pybind11::value_error( "Invalid index: " + std::to_string( i ) + " (must be in [0, " + std::to_string( size ) + "))." );
        return refElement.position( i, c );
      }


      // referenceElementType
      // --------------------

      template< class RefElement >
      inline static GeometryType referenceElementType ( const RefElement &refElement, int i, int c )
      {
        const int size = detail::referenceElementSize( refElement, c );
        if( (i < 0) || (i >= size) )
          throw pybind11::value_error( "Invalid index: " + std::to_string( i ) + " (must be in [0, " + std::to_string( size ) + "))." );
        return refElement.type( i, c );
      }

    } // namespace detail


    template< class RefElement, class... options >
    void registerReferenceElements ( pybind11::module module, pybind11::class_< RefElement, options... > cls )
    {
      using pybind11::operator""_a;

      pybind11::options opts;
      opts.disable_function_signatures();

      static const std::size_t dimension = RefElement::dimension;
      typedef typename RefElement::ctype ctype;
      registerFieldVecMat<FieldVector<ctype,dimension>>::apply();
      registerFieldVecMat<FieldMatrix<ctype,dimension,dimension>>::apply();

      cls.def_property_readonly( "dimension", [] ( pybind11::object ) { return pybind11::int_( RefElement::dimension ); } );

      cls.def( "size", [] ( const RefElement &self, int c ) {
          return detail::referenceElementSize( self, c );
        }, "codim"_a );
      cls.def( "size", [] ( const RefElement &self, int i, int c, int cc ) {
          return detail::referenceElementSize( self, i, c, cc );
         } );
      cls.def( "size", [] ( const RefElement &self, std::tuple< int, int > e, int cc ) {
          return detail::referenceElementSize( self, std::get< 0 >( e ), std::get< 1 >( e ), cc );
        }, "subentity"_a, "codim"_a );

      cls.def( "type", [] ( const RefElement &self, int i, int c ) {
          return detail::referenceElementType( self, i, c );
        }, "index"_a, "codim"_a );
      cls.def( "type", [] ( const RefElement &self, std::tuple< int, int > e ) {
          return detail::referenceElementType( self, std::get< 0 >( e ), std::get< 1 >( e ) );
        }, "subentity"_a );
      cls.def( "types", [] ( const RefElement &self, int c ) {
          const int size = detail::referenceElementSize( self, c );
          pybind11::tuple types( size );
          for( int i = 0; i < size; ++i )
            types[ i ] = pybind11::cast( self.type( i, c ) );
          return types;
        }, "codim"_a,
        R"doc(
          get geometry types of subentities

          Args:
              codim:    codimension of subentities

          Returns:
              tuple containing geometry type for each subentity
        )doc" );

      cls.def( "subEntity", [] ( const RefElement &self, int i, int c, int ii, int cc ) {
          return detail::referenceElementSubEntity( self, i, c, ii, cc );
        } );
      cls.def( "subEntity", [] ( const RefElement &self, std::tuple< int, int > e, std::tuple< int, int > ee ) {
          return detail::referenceElementSubEntity( self, std::get< 0 >( e ), std::get< 1 >( e ), std::get< 0 >( ee ), std::get< 1 >( ee ) );
        } );
      cls.def( "subEntities", [] ( const RefElement &self, int i, int c, int cc ) {
          const int size = detail::referenceElementSize( self, i, c, cc );
          pybind11::tuple subEntities( size );
          for( int ii = 0; ii < size; ++ii )
            subEntities[ ii ] = pybind11::int_( self.subEntity( i, c, ii, cc ) );
          return subEntities;
        } );
      cls.def( "subEntities", [] ( const RefElement &self, std::tuple< int, int > e, int cc ) {
          const int size = detail::referenceElementSize( self, std::get< 0 >( e ), std::get< 1 >( e ), cc );
          pybind11::tuple subEntities( size );
          for( int ii = 0; ii < size; ++ii )
            subEntities[ ii ] = pybind11::int_( self.subEntity( std::get< 0 >( e ), std::get< 1 >( e ), ii, cc ) );
          return subEntities;
        }, "subentity"_a, "codim"_a );

      cls.def( "position", [] ( const RefElement &self, int i, int c ) {
          return detail::referenceElementPosition( self, i, c );
        }, "index"_a, "codim"_a );
      cls.def( "position", [] ( const RefElement &self, std::tuple< int, int > e ) {
          return detail::referenceElementPosition( self, std::get< 0 >( e ), std::get< 1 >( e ) );
        }, "subentity"_a );
      cls.def( "positions", [] ( const RefElement &self, int c ) {
          const int size = detail::referenceElementSize( self, c );
          pybind11::tuple positions( size );
          for( int i = 0; i < size; ++i )
            positions[ i ] = pybind11::cast( self.position( i, c ) );
          return positions;
        }, "codim"_a,
        R"doc(
          get barycenters of subentities

          Args:
              codim:    codimension of subentities

          Returns:
              tuple containing barycenter for each subentity
        )doc" );

#if 0
      // Bug: This property overwrite the method "type"
      cls.def_property_readonly( "type", [] ( const RefElement &self ) {
          return self.type();
        },
        R"doc(
          geometry type of reference element
        )doc" );
#endif
      cls.def_property_readonly( "center", [] ( const RefElement &self ) {
          return self.position( 0, 0 );
        },
        R"doc(
          barycenter of reference domain
        )doc" );
      cls.def_property_readonly( "corners", [] ( const RefElement &self ) {
          const int size = self.size( RefElement::dimension );
          pybind11::tuple corners( size );
          for( int i = 0; i < size; ++i )
            corners[ i ] = pybind11::cast( self.position( i, RefElement::dimension ) );
          return corners;
        },
        R"doc(
          corners of reference domain
        )doc" );
      cls.def_property_readonly( "volume", [] ( const RefElement &self ) {
          return self.volume();
        },
        R"doc(
          volume of reference domain
        )doc" );

      if( RefElement::dimension > 0 )
      {
        cls.def( "integrationOuterNormal", [] ( const RefElement &self, int i ) {
            if( (i < 0) || (i >= self.size( 1 )) )
              throw pybind11::value_error( "Invalid index: " + std::to_string( i ) + " (must be in [0, " + std::to_string( self.size( 1 ) ) + "))." );
            return self.integrationOuterNormal( i );
          }, "index"_a );
        cls.def_property_readonly( "integrationOuterNormals", [] ( const RefElement &self ) {
            const int size = self.size( 1 );
            pybind11::tuple integrationOuterNormals( size );
            for( int i = 0; i < size; ++i )
              integrationOuterNormals[ i ] = pybind11::cast( self.integrationOuterNormal( i ) );
            return integrationOuterNormals;
          } );
      }

      registerQuadratureRule<typename RefElement::ctype, RefElement::dimension>( module );

      module.def( "general", [] ( const GeometryType &gt ) -> pybind11::object {
          if( gt.isNone() )
            return pybind11::none();
          else
            return pybind11::cast( Dune::ReferenceElements< typename RefElement::ctype, RefElement::dimension >::general( gt ), pybind11::return_value_policy::reference );
        } );
      module.def( "rule", [] (const GeometryType &gt, int order) {
            return Dune::QuadratureRules< typename RefElement::ctype, RefElement::dimension >::rule( gt, order );
      }, pybind11::return_value_policy::reference );
    }

    template <int dim>
    void registerReferenceElementToModule( pybind11::module module )
    {
      using pybind11::operator""_a;
      pybind11::module cls0 = module;
      {
        using DuneType = Dune::Geo::ReferenceElement<Dune::Geo::ReferenceElementImplementation<double,dim> >;
        std::string typeName ("Dune::Geo::ReferenceElement<Dune::Geo::ReferenceElementImplementation<double,");
        typeName += std::to_string( dim ) + "> >";
        auto cls = Dune::Python::insertClass< DuneType >( cls0, "ReferenceElements",
              Dune::Python::GenerateTypeName( typeName.c_str() ),
              Dune::Python::IncludeFiles{"dune/python/geometry/referenceelements.hh"}).first;
        Dune::Python::registerReferenceElements( cls0, cls );
      }
    }


  } // namespace Python

} // namespace Dune

#endif // #ifndef DUNE_PYTHON_GEOMETRY_REFERENCEELEMENTS_HH