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
|
// -*- tab-width: 4; indent-tabs-mode: nil; c-basic-offset: 2 -*-
// vi: set et ts=4 sw=2 sts=2:
// SPDX-FileCopyrightInfo: 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_QUADRATURERULES_HH
#define DUNE_PYTHON_GEOMETRY_QUADRATURERULES_HH
#include <array>
#include <tuple>
#include <dune/common/visibility.hh>
#include <dune/geometry/quadraturerules.hh>
#include <dune/geometry/type.hh>
#include <dune/python/common/typeregistry.hh>
#include <dune/python/pybind11/pybind11.h>
#include <dune/python/pybind11/numpy.h>
namespace Dune
{
namespace Python
{
template <class Rule>
auto quadratureToNumpy(const Rule &rule)
{
pybind11::array_t< double > p(
{ static_cast< ssize_t >( Rule::d ), static_cast< ssize_t >( rule.size() ) },
{
static_cast< ssize_t >( sizeof(rule[ 0 ].position()[ 0 ]) ),
static_cast< ssize_t >( sizeof(rule[ 0 ]) )
},
&rule[ 0 ].position()[ 0 ]
);
pybind11::array_t< double > w(
{ static_cast< ssize_t >( rule.size() ) },
{ static_cast< ssize_t >( sizeof(rule[ 0 ]) ) },
&rule[ 0 ].weight()
);
return std::make_pair( p, w );
}
template <class Rule>
auto quadratureToNumpy(pybind11::object self)
{
return quadratureToNumpy( pybind11::cast< const Rule & >( self ) );
}
namespace detail
{
// registerQuadraturePoint
// -----------------------
template< class QP >
inline void registerQuadraturePoint ( pybind11::object scope, pybind11::class_<QP> cls )
{
using pybind11::operator""_a;
typedef typename QP::Vector Vector;
typedef typename QP::Field Field;
cls.def( pybind11::init( [] ( const Vector &x, Field w ) { return new QP( x, w ); } ), "position"_a, "weight"_a );
cls.def_property_readonly( "position", []( const QP &qp ) -> Vector { return qp.position(); } );
cls.def_property_readonly( "weight", &QP::weight );
}
// registerQuadratureRule
// ----------------------
template< class Rule, class... options >
inline void registerQuadratureRule ( pybind11::object scope,
pybind11::class_<Rule,options...> cls )
{
cls.def_property_readonly( "order", &Rule::order );
cls.def_property_readonly( "type", &Rule::type );
cls.def( "get", [] ( pybind11::object self ) {
return quadratureToNumpy<Rule>(self);
} );
cls.def( "__iter__", [] ( const Rule &rule ) { return pybind11::make_iterator( rule.begin(), rule.end() ); } );
}
} // namespace detail
// registerQuadratureRule
// ----------------------
template< class ctype, int dim >
inline auto registerQuadratureRule ( pybind11::object scope )
{
typedef typename Dune::QuadraturePoint< ctype, dim > QP;
auto quadPointCls = insertClass< QP >( scope, "QuadraturePoint",
GenerateTypeName("Dune::QuadratePoint",MetaType<ctype>(),dim),
IncludeFiles{"dune/python/geometry/quadraturerules.hh"});
if (quadPointCls.second)
detail::registerQuadraturePoint( scope, quadPointCls.first );
typedef typename Dune::QuadratureRule< ctype, dim > Rule;
auto quadRule = insertClass< Rule >(scope, "QuadratureRule" + std::to_string(dim),
GenerateTypeName("Dune::QuadrateRule",MetaType<ctype>(),dim),
IncludeFiles{"dune/python/geometry/quadraturerules.hh"});
if (quadRule.second)
detail::registerQuadratureRule( scope, quadRule.first );
return quadRule.first;
}
template< class ctype, int ... dim >
inline auto registerQuadratureRule ( pybind11::object scope, std::integer_sequence< int, dim ... > )
{
return std::make_tuple( registerQuadratureRule< ctype >( scope, std::integral_constant< int, dim >() )... );
}
} // namespace Python
} // namespace Dune
#endif // #ifndef DUNE_PYTHON_GEOMETRY_QUADRATURERULES_HH
|