DNDSR 0.2.1
Distributed Numeric Data Structure for CFV
Loading...
Searching...
No Matches
ArrayEigenMatrixBatch_bind.hpp
Go to the documentation of this file.
1#pragma once
2
3#include <utility>
4
6#include "../Array_bind.hpp"
7
8namespace DNDS
9{
10
12 {
13 return "ArrayEigenMatrixBatch";
14 }
15
17 {
18 return "ArrayEigenMatrixBatchPair";
19 }
20
22
24}
25
26namespace DNDS
27{
28 inline void pybind11_ArrayEigenMatrixBatch_setitem_row(ArrayEigenMatrixBatch &self, index i, const py::list &matList)
29 {
30 using tElem = real;
31 using tReadMap = Eigen::Map<
32 const Eigen::Matrix<tElem, Eigen::Dynamic, Eigen::Dynamic, Eigen::DontAlign | Eigen::ColMajor>,
33 Eigen::Unaligned,
34 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>>;
35 std::vector<tReadMap> mat_maps;
36 for (const auto &v : matList)
37 {
38 if (!py::isinstance<py::buffer>(v))
39 throw std::runtime_error("All elements must be buffer-compatible objects.");
40 auto buf = v.cast<py::buffer>();
41 auto buf_info = buf.request(false);
42 DNDS_assert(buf_info.item_type_is_equivalent_to<tElem>());
43 DNDS_assert_info(buf_info.shape.size() == 2, "need to pass a 2-d array");
44 auto *buf_start_ptr = reinterpret_cast<tElem *>(buf_info.ptr);
45 DNDS_assert(buf_info.strides.size() == 2);
46
47 auto c_mat_map = tReadMap(
48 buf_start_ptr,
49 buf_info.shape[0],
50 buf_info.shape[1],
51 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>(buf_info.strides[1] / sizeof(tElem) /*col stride*/, buf_info.strides[0] / sizeof(tElem) /*row stride*/));
52 mat_maps.push_back(c_mat_map);
53 }
54 self.InitializeWriteRow(i, mat_maps);
55 }
56
58 {
59 using tElem = real;
60
61 auto rowBatch = self[index_];
62 py::list ret;
63 for (int iMat = 0; iMat < rowBatch.Size(); iMat++)
64 {
65 auto mat = rowBatch[iMat];
66 ret.append(
67 py::memoryview::from_buffer<tElem>(
68 mat.data(),
69 {mat.rows(), mat.cols()},
70 {sizeof(tElem) * mat.rowStride(), sizeof(tElem) * mat.colStride()},
71 false)); //! warning: deleting list and Array could make the items dangling
72 }
73 return ret;
74 }
75
76 inline auto pybind11_ArrayEigenMatrixBatch_getitem(ArrayEigenMatrixBatch &self, std::tuple<index, rowsize> index_)
77 {
78 using tElem = real;
79
80 auto mat = self(std::get<0>(index_), std::get<1>(index_));
81 return py::memoryview::from_buffer<tElem>(
82 mat.data(),
83 {mat.rows(), mat.cols()},
84 {sizeof(tElem) * mat.rowStride(), sizeof(tElem) * mat.colStride()},
85 false);
86 }
87
88 inline void pybind11_ArrayEigenMatrixBatch_setitem(ArrayEigenMatrixBatch &self, std::tuple<index, rowsize> index_, const py::buffer &row)
89 {
90 using tElem = real;
91 auto row_info = row.request(false);
92 DNDS_assert(row_info.item_type_is_equivalent_to<tElem>());
93 auto mat = self(std::get<0>(index_), std::get<1>(index_));
94 DNDS_assert_info(row_info.shape.size() == 2, "need to pass a 2-d array");
95 DNDS_assert_info(row_info.shape[0] == mat.rows(), "row size not matching");
96 DNDS_assert_info(row_info.shape[1] == mat.cols(), "col size not matching");
97
98 auto *row_start_ptr = reinterpret_cast<tElem *>(row_info.ptr);
99 DNDS_assert(row_info.strides.size() == 2);
100 auto row_mat_map = Eigen::Map<
101 const Eigen::Matrix<tElem, Eigen::Dynamic, Eigen::Dynamic, Eigen::DontAlign | Eigen::ColMajor>,
102 Eigen::Unaligned,
103 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>>(
104 row_start_ptr,
105 row_info.shape[0],
106 row_info.shape[1],
107 Eigen::Stride<Eigen::Dynamic, Eigen::Dynamic>(row_info.strides[1] / sizeof(tElem) /*col stride*/, row_info.strides[0] / sizeof(tElem) /*row stride*/));
108 mat = row_mat_map;
109 }
110}
111
112namespace DNDS
113{
116 {
117 auto ParArray_ = pybind11_ParArray_get_class<real, NonUniformSize>(m);
118 return {m, pybind11_ArrayEigenMatrixBatch_name().c_str(), ParArray_};
119 }
120
123 {
124 return {m.attr(pybind11_ArrayEigenMatrixBatch_name().c_str())};
125 }
126
127 inline void pybind11_ArrayEigenMatrixBatch_define(py::module_ &m)
128 {
129
130 using TArrayEigenMatrixBatch = ArrayEigenMatrixBatch;
131 auto ArrayEigenMatrixBatch_ = pybind11_ArrayEigenMatrixBatch_declare(m);
132
133 ArrayEigenMatrixBatch_
134 // we only bind the non-default ctor here
135 .def(py::init<const MPIInfo &>(), py::arg("nmpi"))
136 .def(
137 "BatchSize",
138 [](TArrayEigenMatrixBatch &self, index i)
139 {
140 return self.BatchSize(i);
141 },
142 py::arg("iRow"))
143 .def(
144 "InitializeWriteRow",
145 [](TArrayEigenMatrixBatch &self, index i, const py::list &matList)
146 {
148 },
149 py::arg("i"), py::arg("matList"));
150 ArrayEigenMatrixBatch_
151 .def("clone", [](TArrayEigenMatrixBatch &self)
152 {
153 auto arr = std::make_shared<TArrayEigenMatrixBatch>(self);
154 return arr; });
155 ArrayEigenMatrixBatch_
156 .def(
157 "__getitem__",
158 [](TArrayEigenMatrixBatch &self, index index_)
159 {
161 },
162 py::keep_alive<0, 1>())
163 .def(
164 "__getitem__",
165 [](TArrayEigenMatrixBatch &self, std::tuple<index, rowsize> index_)
166 {
167 return pybind11_ArrayEigenMatrixBatch_getitem(self, index_);
168 },
169 py::keep_alive<0, 1>())
170 .def(
171 "__setitem__",
172 [](TArrayEigenMatrixBatch &self, std::tuple<index, rowsize> index_, const py::buffer &row)
173 {
175 });
176
177 ArrayEigenMatrixBatch_
178 .def("to_device", [](TArrayEigenMatrixBatch &self, const std::string &backend)
179 { self.to_device(device_backend_name_to_enum(backend)); }, py::arg("backend"))
180 .def("to_host", &TArrayEigenMatrixBatch::to_host);
181 }
182}
183
184namespace DNDS
185{
191
194 {
195 return {m.attr(pybind11_ArrayEigenMatrixBatchPair_name().c_str())};
196 }
197
199 {
200
201 using TPair = ArrayEigenMatrixBatchPair;
203
204 pybind11_ArrayPairGenericBindBasics<TPair>(Pair_);
205
206 Pair_
207 .def("RowSize", [](const TPair &self, index i)
208 { return self.RowSize(i); }, py::arg("i"));
209 Pair_
210 .def(
211 "__getitem__",
212 [](TPair &self, index index_)
213 {
214 return self.runFunctionAppendedIndex(index_, [&](auto &ar, index iC) //*note the auto&& reference here!!!
216 },
217 py::keep_alive<0, 1>())
218 .def(
219 "InitializeWriteRow",
220 [](TPair &self, index index_, const py::buffer &row)
221 {
222 self.runFunctionAppendedIndex(index_, [&](auto &ar, index iC) //*note the auto&& reference here!!!
224 });
225
226 Pair_
227 .def(
228 "__getitem__",
229 [](TPair &self, std::tuple<index, rowsize> index_)
230 {
231 return self.runFunctionAppendedIndex(std::get<0>(index_), [&](auto &ar, index iC) //*note the auto&& reference here!!!
232 { return pybind11_ArrayEigenMatrixBatch_getitem(ar, std::make_tuple(iC, std::get<1>(index_))); });
233 },
234 py::keep_alive<0, 1>())
235 .def(
236 "__setitem__",
237 [](TPair &self, std::tuple<index, rowsize> index_, const py::buffer &row)
238 {
239 self.runFunctionAppendedIndex(std::get<0>(index_), [&](auto &ar, index iC) //*note the auto&& reference here!!!
240 { pybind11_ArrayEigenMatrixBatch_setitem(ar, std::make_tuple(iC, std::get<1>(index_)), row); });
241 });
242 }
243}
244
245namespace DNDS
246{
247
253}
Batch of variable-sized Eigen matrices stored in CSR layout.
#define DNDS_assert_info(expr, info)
Debug-only assertion with an extra std::string info message.
Definition Errors.hpp:117
#define DNDS_assert(expr)
Debug-only assertion (compiled out when DNDS_NDEBUG is defined). Prints the expression + file/line + ...
Definition Errors.hpp:112
Eigen::Matrix< real, 3, 3 > m
CSR array storing a variable-sized batch of Eigen matrices per row.
void InitializeWriteRow(index i, const std::vector< t_matrices_elem > &matrices)
the host side operators are provided as implemented
void pybind11_ArrayEigenMatrixBatchPair_define(py::module_ &m)
void pybind11_ArrayEigenMatrixBatch_setitem_row(ArrayEigenMatrixBatch &self, index i, const py::list &matList)
py_class_ssp< ArrayEigenMatrixBatchPair > tPy_ArrayEigenMatrixBatchPair
py::classh< T > py_class_ssp
tPy_ArrayEigenMatrixBatchPair pybind11_ArrayEigenMatrixBatchPair_declare(py::module_ &m)
std::string pybind11_ArrayEigenMatrixBatch_name()
DeviceBackend device_backend_name_to_enum(std::string_view s)
Inverse of device_backend_name. Returns Unknown for unrecognised names.
tPy_ArrayEigenMatrixBatch pybind11_ArrayEigenMatrixBatch_declare(py::module_ &m)
py_class_ssp< ArrayEigenMatrixBatch > tPy_ArrayEigenMatrixBatch
tPy_ArrayEigenMatrixBatch pybind11_ArrayEigenMatrixBatch_get_class(py::module_ &m)
void pybind11_ArrayEigenMatrixBatch_define(py::module_ &m)
auto pybind11_ArrayEigenMatrixBatch_getitem(ArrayEigenMatrixBatch &self, std::tuple< index, rowsize > index_)
int64_t index
Global row / DOF index type (signed 64-bit; handles multi-billion-cell meshes).
Definition Defines.hpp:112
ArrayPair< ArrayEigenMatrixBatch > ArrayEigenMatrixBatchPair
ArrayPair alias for per-row variable-size Eigen matrix batches.
double real
Canonical floating-point scalar used throughout DNDSR (double precision).
Definition Defines.hpp:110
tPy_ArrayEigenMatrixBatchPair pybind11_ArrayEigenMatrixBatchPair_get_class(py::module_ &m)
void pybind11_ArrayEigenMatrixBatch_setitem(ArrayEigenMatrixBatch &self, std::tuple< index, rowsize > index_, const py::buffer &row)
std::string pybind11_ArrayEigenMatrixBatchPair_name()
auto pybind11_ArrayEigenMatrixBatch_getitem_row(ArrayEigenMatrixBatch &self, index index_)
void pybind11_bind_ArrayEigenMatrixBatch_All(py::module_ &m)
Eigen::Matrix< real, 5, 1 > v