17#include <vtkDataArray.h>
27class vtkUnstructuredGrid;
44 vtkUnstructuredGrid* grid,
45 bool const compute_element_neighbors =
false,
46 std::string
const& mesh_name =
"vtkUnstructuredGrid");
53 std::vector<double>
const& elevation,
54 std::vector<int>& node_idx_map,
60 std::vector<double>
const& pix_val,
61 std::vector<bool>
const& pix_vis,
62 std::vector<MeshLib::Node*>
const& nodes,
63 std::vector<int>
const& node_idx_map,
64 std::size_t
const imgHeight,
65 std::size_t
const imgWidth,
71 std::vector<double>
const& pix_val,
72 std::vector<bool>
const& pix_vis,
73 const std::size_t& imgHeight,
74 const std::size_t& imgWidth,
77 for (std::size_t i = 0; i < imgHeight; i++)
79 for (std::size_t j = 0; j < imgWidth; j++)
81 if (!pix_vis[i * imgWidth + j])
85 T val(
static_cast<T
>(pix_val[i * (imgWidth + 1) + j]));
88 prop_vec.push_back(val);
89 prop_vec.push_back(val);
93 prop_vec.push_back(val);
103 template <
typename T>
116 if (
sizeof(T) != array.GetDataTypeSize())
119 "Trying to convert a vtk array '{:s}' with data type '{:s}' of "
120 "size {:d} to a different sized type '{:s}' of size {:d}.",
121 array.GetName(), array.GetDataTypeAsString(),
122 array.GetDataTypeSize(),
typeid(T).name(),
sizeof(T));
125 vtkIdType
const nTuples(array.GetNumberOfTuples());
126 int const nComponents(array.GetNumberOfComponents());
127 char const*
const array_name(array.GetName());
130 array_name, type, nComponents);
133 WARN(
"Array {:s} could not be converted to PropertyVector.",
137 vec->reserve(nTuples * nComponents);
138 auto* data_array =
static_cast<T*
>(array.GetVoidPointer(0));
139 std::copy(&data_array[0],
140 &data_array[nTuples * nComponents],
141 std::back_inserter(*vec));
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition of mesh-related Enumerations.
Definition of the class Properties that implements a container of properties.
Definition of the GeoLib::Raster class.
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
PropertyVector< T > * createNewPropertyVector(std::string_view name, MeshItemType mesh_item_type, std::size_t n_components=1)
Converter for VtkUnstructured Grids to OGS meshes.
static std::vector< MeshLib::Element * > createElementVector(std::vector< double > const &pix_val, std::vector< bool > const &pix_vis, std::vector< MeshLib::Node * > const &nodes, std::vector< int > const &node_idx_map, std::size_t const imgHeight, std::size_t const imgWidth, MeshElemType elem_type)
Creates a mesh element vector based on image data.
static void fillPropertyVector(MeshLib::PropertyVector< T > &prop_vec, std::vector< double > const &pix_val, std::vector< bool > const &pix_vis, const std::size_t &imgHeight, const std::size_t &imgWidth, MeshElemType elem_type)
Creates a scalar array/mesh property based on pixel values.
static void convertTypedArray(vtkDataArray &array, MeshLib::Properties &properties, MeshLib::MeshItemType type)
static void convertScalarArrays(vtkUnstructuredGrid &grid, MeshLib::Mesh &mesh)
static MeshLib::Mesh * convertUnstructuredGrid(vtkUnstructuredGrid *grid, bool const compute_element_neighbors=false, std::string const &mesh_name="vtkUnstructuredGrid")
Converts a vtkUnstructuredGrid object to a Mesh.
static void convertArray(vtkDataArray &array, MeshLib::Properties &properties, MeshLib::MeshItemType type)
static std::vector< MeshLib::Node * > createNodeVector(std::vector< double > const &elevation, std::vector< int > &node_idx_map, GeoLib::RasterHeader const &header, bool use_elevation)
MeshElemType
Types of mesh elements supported by OpenGeoSys. Values are from VTKCellType enum.