File: densevector.hh

package info (click to toggle)
dune-common 2.10.0-6
  • links: PTS, VCS
  • area: main
  • in suites: forky, sid
  • size: 5,824 kB
  • sloc: cpp: 52,256; python: 3,979; sh: 1,658; makefile: 17
file content (194 lines) | stat: -rw-r--r-- 7,933 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
// -*- tab-width: 2; 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_COMMON_DENSEVECTOR_HH
#define DUNE_PYTHON_COMMON_DENSEVECTOR_HH

#include <string>
#include <tuple>
#include <type_traits>
#include <utility>

#include <dune/python/common/vector.hh>
#include <dune/python/pybind11/pybind11.h>
#include <dune/python/pybind11/operators.h>

namespace Dune
{

  namespace Python
  {

    namespace detail
    {

      // registerScalarCopyingDenseVectorMethods
      // ---------------------------------------

      template< class T, class... options >
      inline static std::enable_if_t< std::is_copy_constructible< T >::value && (T::dimension == 1) >
      registerScalarCopyingDenseVectorMethods ( pybind11::class_< T, options... > cls, PriorityTag< 2 > )
      {
        using ValueType = typename T::value_type;

        cls.def( "__add__", [] ( const T &self, int a ) { T *copy = new T( self ); (*copy)[ 0 ] += ValueType( a ); return copy; } );
        cls.def( "__add__", [] ( const T &self, const ValueType &a ) { T *copy = new T( self ); (*copy)[ 0 ] += a; return copy; } );
        cls.def( "__sub__", [] ( const T &self, int a ) { T *copy = new T( self ); (*copy)[ 0 ] -= ValueType( a ); return copy; } );
        cls.def( "__sub__", [] ( const T &self, const ValueType &a ) { T *copy = new T( self ); (*copy)[ 0 ] -= a; return copy; } );

        cls.def( "__radd__", [] ( const T &self, int a ) { T *copy = new T( self ); (*copy)[ 0 ] = ValueType( a ) + (*copy)[ 0 ]; return copy; } );
        cls.def( "__radd__", [] ( const T &self, const ValueType &a ) { T *copy = new T( self ); (*copy)[ 0 ] = a + (*copy)[ 0 ]; return copy; } );
        cls.def( "__rsub__", [] ( const T &self, int a ) { T *copy = new T( self ); (*copy)[ 0 ] = ValueType( a ) - (*copy)[ 0 ]; return copy; } );
        cls.def( "__rsub__", [] ( const T &self, const ValueType &a ) { T *copy = new T( self ); (*copy)[ 0 ] = a - (*copy)[ 0 ]; return copy; } );
      }

      template< class T, class... options >
      inline static std::enable_if_t< std::is_copy_constructible< T >::value >
      registerScalarCopyingDenseVectorMethods ( pybind11::class_< T, options... > cls, PriorityTag< 1 > )
      {
        using ValueType = typename T::value_type;

        cls.def( "__add__", [] ( pybind11::object self, int a ) {
            if( a != 0 )
              throw pybind11::value_error( "Cannot add " + std::to_string( a ) + " to multidimensional dense vector." );
            return self;
          } );
        cls.def( "__sub__", [] ( pybind11::object self, int a ) {
            if( a != 0 )
              throw pybind11::value_error( "Cannot subtract " + std::to_string( a ) + " from multidimensional dense vector." );
            return self;
          } );

        cls.def( "__radd__", [] ( pybind11::object self, int a ) {
            if( a != 0 )
              throw pybind11::value_error( "Cannot add multidimensional dense vector to " + std::to_string( a ) + "." );
            return self;
          } );
        cls.def( "__rsub__", [] ( const T &self, int a ) {
            if( a != 0 )
              throw pybind11::value_error( "Cannot subtract multidimensional dense vector from " + std::to_string( a ) + "." );
            T *copy = new T( self ); *copy *= ValueType( -1 ); return copy;
          } );
      }

      template< class T, class... options >
      inline static void registerScalarCopyingDenseVectorMethods ( pybind11::class_< T, options... > cls, PriorityTag< 0 > )
      {}

      template< class T, class... options >
      inline static void registerScalarCopyingDenseVectorMethods ( pybind11::class_< T, options... > cls )
      {
        registerScalarCopyingDenseVectorMethods ( cls, PriorityTag< 42 >() );
      }




      // registerCopyingDenseVectorMethods
      // ---------------------------------

      template< class T, class... options >
      inline static std::enable_if_t< std::is_copy_constructible< T >::value >
      registerCopyingDenseVectorMethods ( pybind11::class_< T, options... > cls, PriorityTag< 1 > )
      {
        using ValueType = typename T::value_type;

        using pybind11::operator""_a;

        cls.def( "__pos__", [] ( pybind11::object self ) { return self; } );
        cls.def( "__neg__", [] ( T &self ) { T *copy = new T( self ); *copy *= ValueType( -1 ); return copy; } );

        cls.def( pybind11::self + pybind11::self );
        cls.def( pybind11::self - pybind11::self );

        cls.def( "__add__", [] ( T &self, pybind11::list x ) { return self + x.cast< T >(); }, "x"_a );
        cls.def( "__sub__", [] ( T &self, pybind11::list x ) { return self - x.cast< T >(); }, "x"_a );

        cls.def( "__radd__", [] ( T &self, pybind11::list x ) { return x.cast< T >() + self; }, "x"_a );
        cls.def( "__rsub__", [] ( T &self, pybind11::list x ) { return x.cast< T >() - self; }, "x"_a );

        cls.def( "__mul__", [] ( const T &self, ValueType x ) { T *copy = new T( self ); *copy *= x; return copy; }, "x"_a );
        cls.def( "__div__", [] ( const T &self, ValueType x ) { T *copy = new T( self ); *copy /= x; return copy; }, "x"_a );
        cls.def( "__truediv__", [] ( const T &self, ValueType x ) { T *copy = new T( self ); *copy /= x; return copy; }, "x"_a );

        cls.def( "__rmul__", [] ( const T &self, ValueType x ) { T *copy = new T( self ); *copy *= x; return copy; }, "x"_a );
      }

      template< class T, class... options >
      inline static void registerCopyingDenseVectorMethods ( pybind11::class_< T, options... > cls, PriorityTag< 0 > )
      {}

      template< class T, class... options >
      inline static void registerCopyingDenseVectorMethods ( pybind11::class_< T, options... > cls )
      {
        registerCopyingDenseVectorMethods ( cls, PriorityTag< 42 >() );
      }

    } // namespace detail



    // registerDenseVector
    // -------------------

    template< class T, class... options >
    inline static void registerDenseVector ( pybind11::class_< T, options... > cls )
    {
      using ValueType = typename T::value_type;

      using pybind11::operator""_a;

      cls.def( "assign", [] ( T &self, const T &x ) { self = x; }, "x"_a );

      cls.def( "__getitem__", [] ( const T &self, std::size_t i ) -> ValueType {
          if( i < self.size() )
            return self[ i ];
          else
            throw pybind11::index_error();
        }, "i"_a );

      cls.def( "__setitem__", [] ( T &self, std::size_t i, ValueType x ) {
          if( i < self.size() )
            self[ i ] = x;
          else
            throw pybind11::index_error();
        }, "i"_a, "x"_a );

      cls.def( "__len__", [] ( const T &self ) -> std::size_t { return self.size(); } );

      cls.def( pybind11::self += pybind11::self );

// silence a warning (false positive) emitted by clang
// https://bugs.llvm.org/show_bug.cgi?id=43124
#ifdef __clang__
#pragma GCC diagnostic push
#pragma GCC diagnostic ignored "-Wself-assign-overloaded"
#endif

      cls.def( pybind11::self -= pybind11::self );

#ifdef __clang__
#pragma GCC diagnostic pop
#endif

      cls.def( pybind11::self == pybind11::self );
      cls.def( pybind11::self != pybind11::self );

      cls.def( pybind11::self += ValueType() );
      cls.def( pybind11::self -= ValueType() );
      cls.def( pybind11::self *= ValueType() );
      cls.def( pybind11::self /= ValueType() );

      detail::registerOneTensorInterface( cls );
      detail::registerCopyingDenseVectorMethods( cls );
      detail::registerScalarCopyingDenseVectorMethods( cls );

      pybind11::implicitly_convertible< pybind11::list, T >();
    }

  } // namespace Python

} // namespace Dune

#endif // #ifndef DUNE_PYTHON_COMMON_DENSEVECTOR_HH