OGS
OGSMesh.h
Go to the documentation of this file.
1// SPDX-FileCopyrightText: Copyright (c) OpenGeoSys Community (opengeosys.org)
2// SPDX-License-Identifier: BSD-3-Clause
3
4#include <pybind11/numpy.h>
5#include <pybind11/pybind11.h>
6
7#include "MeshLib/Mesh.h"
8#include "MeshLib/MeshEnums.h"
9
10// Needs to be exported, see
11// https://pybind11.readthedocs.io/en/stable/advanced/misc.html#partitioning-code-over-multiple-extension-modules
12class PYBIND11_EXPORT OGSMesh
13{
14public:
15 explicit OGSMesh(MeshLib::Mesh& mesh);
16
17 std::vector<double> getPointCoordinates() const;
18 std::pair<std::vector<int>, std::vector<int>> getCells() const;
19
20 std::vector<std::string> getDataArrayNames() const;
21
22 MeshLib::MeshItemType meshItemType(std::string_view const name) const;
23
24 template <typename T>
25 pybind11::array& dataArray(std::string const& name)
26 {
27 auto const data_array_it = data_array_mapping.find(name);
28 if (data_array_it != data_array_mapping.end())
29 {
30 INFO("found data array '{}' with address: {}", name,
31 fmt::ptr(&(data_array_it->second)));
32 return data_array_it->second;
33 }
34
35 auto& mesh_properties = _mesh.getProperties();
36
37 // check if PropertyVector with specified data type T exists
38 auto* property = mesh_properties.getPropertyVector<T>(name);
39 if (property == nullptr)
40 {
41 OGS_FATAL("Couldn't access data array '{}'.", name);
42 }
43
44 // Capsule ties lifetime of the numpy array to the PropertyVector
45 auto capsule = pybind11::capsule(property, [](void* /*ignored*/) {});
46
47 auto const n_components = property->getNumberOfGlobalComponents();
48 if (n_components == 1)
49 {
50 auto const& [it, success] = data_array_mapping.insert(
51 {name,
52 pybind11::array(
53 pybind11::buffer_info(
54 property->data(), // data pointer
55 sizeof(T), // item size
56 pybind11::format_descriptor<T>::format(), // dtype
57 // format
58 // string
59 1, // ndim
60 {property->size()}, // shape
61 {sizeof(T)} // strides
62 ),
63 capsule)});
64 if (!success)
65 {
66 OGS_FATAL("Could not insert data array '{}' into internal map.",
67 name);
68 }
69 INFO("insert data array '{}' with address: {}", name,
70 fmt::ptr(&(it->second)));
71 return it->second;
72 }
73
74 // 2D case: (n_items, n_components)
75 auto const& [it, success] = data_array_mapping.insert(
76 {name,
77 pybind11::array(
78 pybind11::buffer_info(
79 property->data(), // data pointer
80 sizeof(T), // item size
81 pybind11::format_descriptor<T>::format(), // dtype
82 // format
83 // string
84 2, // ndim
85 std::vector<pybind11::ssize_t>{
86 static_cast<pybind11::ssize_t>(property->size() /
87 n_components),
88 static_cast<pybind11::ssize_t>(
89 n_components)}, // shape
90 std::vector<pybind11::ssize_t>{
91 static_cast<pybind11::ssize_t>(sizeof(T) *
92 n_components),
93 static_cast<pybind11::ssize_t>(sizeof(T))} // strides
94 ),
95 capsule)});
96 if (!success)
97 {
98 OGS_FATAL("Could not insert data array '{}' into internal map.",
99 name);
100 }
101 INFO("insert data array '{}' with address: {}", name,
102 fmt::ptr(&(it->second)));
103 return it->second;
104 }
105
106 pybind11::object dataArray_dispatch(std::string const& name,
107 std::string const& dtype)
108 {
109 if (dtype == "double")
110 {
111 return dataArray<double>(name);
112 }
113 if (dtype == "float")
114 {
115 return dataArray<float>(name);
116 }
117 if (dtype == "int")
118 {
119 return dataArray<int>(name);
120 }
121 if (dtype == "int64")
122 {
123 return dataArray<int64_t>(name);
124 }
125 if (dtype == "std::size_t")
126 {
127 return dataArray<std::size_t>(name);
128 }
129 if (dtype == "char")
130 {
131 return dataArray<char>(name);
132 }
133
134 throw std::runtime_error("Unsupported dtype: " + dtype);
135 }
136
137private:
139 std::map<std::string, pybind11::array> data_array_mapping;
140};
#define OGS_FATAL(...)
Definition Error.h:19
void INFO(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:28
MeshLib::MeshItemType meshItemType(std::string_view const name) const
Definition OGSMesh.cpp:48
pybind11::array & dataArray(std::string const &name)
Definition OGSMesh.h:25
MeshLib::Mesh & _mesh
Definition OGSMesh.h:138
std::vector< double > getPointCoordinates() const
Definition OGSMesh.cpp:23
pybind11::object dataArray_dispatch(std::string const &name, std::string const &dtype)
Definition OGSMesh.h:106
std::vector< std::string > getDataArrayNames() const
Definition OGSMesh.cpp:43
std::map< std::string, pybind11::array > data_array_mapping
Definition OGSMesh.h:139
OGSMesh(MeshLib::Mesh &mesh)
Definition OGSMesh.cpp:21
std::pair< std::vector< int >, std::vector< int > > getCells() const
Definition OGSMesh.cpp:29