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
|
// 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_UTILITY_VECTORCOMMDATAHANDLE_HH
#define DUNE_PYTHON_UTILITY_VECTORCOMMDATAHANDLE_HH
#include <cassert>
#include <cstddef>
#include <algorithm>
#include <type_traits>
#include <utility>
#include <dune/common/visibility.hh>
#include <dune/geometry/type.hh>
#include <dune/grid/common/datahandleif.hh>
#include <dune/python/pybind11/numpy.h>
#include <dune/python/pybind11/pybind11.h>
namespace Dune
{
// External Forward Declarations
// -----------------------------
template< class >
class MultipleCodimMultipleGeomTypeMapper;
namespace Python
{
// NumPyCommDataHandle
// -------------------
template< class Mapper, class T, class Function >
class NumPyCommDataHandle;
template< class GV, class T, class Function >
class DUNE_PRIVATE NumPyCommDataHandle< MultipleCodimMultipleGeomTypeMapper< GV >, T, Function >
: public CommDataHandleIF< NumPyCommDataHandle< MultipleCodimMultipleGeomTypeMapper< GV >, T, Function >, T >
{
typedef NumPyCommDataHandle< MultipleCodimMultipleGeomTypeMapper< GV >, T, Function > This;
typedef MultipleCodimMultipleGeomTypeMapper< GV > Mapper;
public:
NumPyCommDataHandle ( const Mapper &mapper, std::vector< pybind11::array_t< T > > arrays, Function function = Function() )
: mapper_( mapper ), buffers_( arrays.size() ), function_( function )
{
std::transform( arrays.begin(), arrays.end(), buffers_.begin(), [] ( pybind11::array_t< T > &a ) { return a.request(); } );
itemSize_ = 0;
for( const pybind11::buffer_info &buffer : buffers_ )
{
if( static_cast< std::size_t >( buffer.shape[ 0 ] ) != mapper_.size() )
pybind11::value_error( "Array does not match mapper in construction of NumPyCommDataHandle." );
itemSize_ += std::accumulate( buffer.shape.begin()+1, buffer.shape.end(), std::size_t( 1 ), std::multiplies< std::size_t >() );
}
}
NumPyCommDataHandle ( const Mapper &mapper, pybind11::array_t< T > array, Function function = Function() )
: NumPyCommDataHandle( mapper, std::vector< pybind11::array_t< T > >{ array }, function )
{}
bool contains ( int dim, int codim ) const
{
const auto &types = mapper_.types( codim );
return std::any_of( types.begin(), types.end(), [ this ] ( GeometryType gt ) { return mapper_.size( gt ) > 0; } );
}
bool fixedSize ( int dim, int codim ) const
{
const auto &types = mapper_.types( codim );
return (std::adjacent_find( types.begin(), types.end(), [ this ] ( GeometryType a, GeometryType b ) { return mapper_.size( a ) != mapper_.size( b ); } ) == types.end());
}
template< class Entity >
std::size_t size ( const Entity &entity ) const
{
return mapper_.size( entity.type() ) * itemSize_;
}
template< class CommBuffer, class Entity >
void gather ( CommBuffer &commBuffer, const Entity &entity ) const
{
for( const pybind11::buffer_info &buffer : buffers_ )
for( const auto index : mapper_.indices( entity ) )
gather( commBuffer, buffer, 1, index*buffer.strides[ 0 ] );
}
template< class CommBuffer, class Entity >
void scatter ( CommBuffer &commBuffer, const Entity &entity, std::size_t n )
{
assert( n == size( entity ) );
for( const pybind11::buffer_info &buffer : buffers_ )
for( const auto index : mapper_.indices( entity ) )
scatter( commBuffer, buffer, 1, index*buffer.strides[ 0 ] );
}
private:
template< class CommBuffer >
void gather ( CommBuffer &commBuffer, const pybind11::buffer_info &buffer, ssize_t dim, ssize_t pos ) const
{
if( dim < buffer.ndim )
{
for( ssize_t i = 0; i < buffer.shape[ dim ]; ++i )
gather( commBuffer, buffer, dim+1, pos + i*buffer.strides[ dim ] );
}
else
commBuffer.write( *reinterpret_cast< T * >( static_cast< char * >( buffer.ptr ) + pos ) );
}
template< class CommBuffer >
void scatter ( CommBuffer &commBuffer, const pybind11::buffer_info &buffer, ssize_t dim, ssize_t pos )
{
if( dim < buffer.ndim )
{
for( ssize_t i = 0; i < buffer.shape[ dim ]; ++i )
scatter( commBuffer, buffer, dim+1, pos + i*buffer.strides[ dim ] );
}
else
{
T &local = *reinterpret_cast< T * >( static_cast< char * >( buffer.ptr ) + pos );
T remote;
commBuffer.read( remote );
local = function_( local, remote );
}
}
const Mapper &mapper_;
std::vector< pybind11::buffer_info > buffers_;
std::size_t itemSize_;
Function function_;
};
// vectorCommDataHandle
// --------------------
template< class Mapper, class T, class Function >
inline static NumPyCommDataHandle< Mapper, T, Function >
numPyCommDataHandle ( const Mapper &mapper, pybind11::array_t< T > array, Function function )
{
return NumPyCommDataHandle< Mapper, T, Function >( mapper, std::move( array ), std::move( function ) );
}
template< class Mapper, class T, class Function >
inline static NumPyCommDataHandle< Mapper, T, Function >
numPyCommDataHandle( const Mapper &mapper, std::vector<pybind11::array_t< T >> array, Function function )
{
return NumPyCommDataHandle< Mapper, T, Function >( mapper, std::move( array ), std::move( function ) );
}
template< class Mapper, class T, class Function, class... options >
void registerDataHandle ( pybind11::handle module, pybind11::class_< NumPyCommDataHandle< Mapper, T, Function >, options... > cls )
{
cls.def( pybind11::init( [] ( Mapper &mapper, pybind11::array_t< T > array, Function function ) {
return NumPyCommDataHandle< Mapper, T, Function >( mapper, array, function );
} ), pybind11::keep_alive< 1, 2 >() );
cls.def( pybind11::init( [] ( Mapper &mapper, std::vector< pybind11::array_t< T > > arrays, Function function ) {
return NumPyCommDataHandle< Mapper, T, Function >( mapper, arrays, function );
} ), pybind11::keep_alive< 1, 2 >() );
}
} // namespace Python
} // namespace Dune
#endif // #ifndef DUNE_PYTHON_UTILITY_VECTORCOMMDATAHANDLE_HH
|