File: densematrix.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 (89 lines) | stat: -rw-r--r-- 3,139 bytes parent folder | download | duplicates (2)
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
// 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_DENSEMATRIX_HH
#define DUNE_PYTHON_COMMON_DENSEMATRIX_HH

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

#include <dune/python/pybind11/extensions.h>
#include <dune/python/pybind11/operators.h>
#include <dune/python/pybind11/pybind11.h>

namespace Dune
{

  namespace Python
  {

    // registerDenseMatrix
    // -------------------

    template< class Matrix >
    void registerDenseMatrix ( pybind11::class_< Matrix > cls )
    {
      typedef typename Matrix::field_type field_type;
      typedef typename Matrix::row_type row_type;
      typedef typename Matrix::row_reference row_reference;

      cls.def( "__getitem__", [] ( Matrix &self, std::size_t i ) -> row_reference {
          if( i < self.mat_rows() )
            return self[ i ];
          else
            throw pybind11::index_error();
        }, (std::is_reference< row_reference >::value ? pybind11::return_value_policy::reference : pybind11::return_value_policy::move), pybind11::keep_alive< 0, 1 >() );

      cls.def( "__setitem__", [] ( Matrix &self, std::size_t i, pybind11::object l ) {
          if( i < self.mat_rows() )
          {
            row_type v = l.cast< row_type >();
            std::size_t size = std::min( self.mat_cols(), v.size() );

            for( std::size_t j = 0; j < size; ++j )
              self[ i ][ j ] = v[ j ];
          }
          else
            throw pybind11::index_error();
        } );

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

      cls.def( "invert", [] ( Matrix &self ) { self.invert(); } );

      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 *= field_type() );
      cls.def( pybind11::self /= field_type() );

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

      cls.def_property_readonly( "frobenius_norm", [] ( const Matrix &self ) { return self.frobenius_norm(); } );
      cls.def_property_readonly( "frobenius_norm2", [] ( const Matrix &self ) { return self.frobenius_norm2(); } );
      cls.def_property_readonly( "infinity_norm", [] ( const Matrix &self ) { return self.infinity_norm(); } );
      cls.def_property_readonly( "infinity_norm_real", [] ( const Matrix &self ) { return self.infinity_norm_real(); } );

      cls.def_property_readonly( "rows", [] ( const Matrix &self ) { return self.mat_rows(); } );
      cls.def_property_readonly( "cols", [] ( const Matrix &self ) { return self.mat_cols(); } );
    }

  } // namespace Python

} // namespace Dune

#endif // #ifndef DUNE_PYTHON_COMMON_DENSEMATRIX_HH