OGS
VtkMeshConverter.h
Go to the documentation of this file.
1
15#pragma once
16
17#include <vtkDataArray.h>
18#include <vtkType.h>
19
20#include "BaseLib/Logging.h"
21#include "GeoLib/Raster.h"
22#include "MeshLib/MeshEnums.h"
23#include "MeshLib/Properties.h"
25
26class vtkUnstructuredGrid; // For conversion vom vtk to ogs mesh
27class vtkDataArray; // For node/cell properties
28
29namespace MeshLib
30{
31class Mesh;
32class Element;
33class Node;
34
39{
40public:
43 vtkUnstructuredGrid* grid,
44 bool const compute_element_neighbors = false,
45 std::string const& mesh_name = "vtkUnstructuredGrid");
46
47private:
48 static void convertScalarArrays(vtkUnstructuredGrid& grid,
49 MeshLib::Mesh& mesh);
50
51 static std::vector<MeshLib::Node*> createNodeVector(
52 std::vector<double> const& elevation,
53 std::vector<int>& node_idx_map,
54 GeoLib::RasterHeader const& header,
55 bool use_elevation);
56
58 static std::vector<MeshLib::Element*> createElementVector(
59 std::vector<double> const& pix_val,
60 std::vector<bool> const& pix_vis,
61 std::vector<MeshLib::Node*> const& nodes,
62 std::vector<int> const& node_idx_map,
63 std::size_t const imgHeight,
64 std::size_t const imgWidth,
65 MeshElemType elem_type);
66
68 template <typename T>
70 std::vector<double> const& pix_val,
71 std::vector<bool> const& pix_vis,
72 const std::size_t& imgHeight,
73 const std::size_t& imgWidth,
74 MeshElemType elem_type)
75 {
76 for (std::size_t i = 0; i < imgHeight; i++)
77 {
78 for (std::size_t j = 0; j < imgWidth; j++)
79 {
80 if (!pix_vis[i * imgWidth + j])
81 {
82 continue;
83 }
84 T val(static_cast<T>(pix_val[i * (imgWidth + 1) + j]));
85 if (elem_type == MeshElemType::TRIANGLE)
86 {
87 prop_vec.push_back(val);
88 prop_vec.push_back(val);
89 }
90 else if (elem_type == MeshElemType::QUAD)
91 {
92 prop_vec.push_back(val);
93 }
94 }
95 }
96 }
97
98 static void convertArray(vtkDataArray& array,
99 MeshLib::Properties& properties,
101
102 template <typename T>
103 static void convertTypedArray(vtkDataArray& array,
104 MeshLib::Properties& properties,
106 {
107 // Keep for debugging purposes, since the vtk array type and the end
108 // type T are not always the same.
109 // DBUG(
110 // "Converting a vtk array '{:s}' with data type '{:s}' of "
111 // "size {:d} to a type '{:s}' of size {:d}.",
112 // array.GetName(), array.GetDataTypeAsString(),
113 // array.GetDataTypeSize(), typeid(T).name(), sizeof(T));
114
115 if (sizeof(T) != array.GetDataTypeSize())
116 {
117 OGS_FATAL(
118 "Trying to convert a vtk array '{:s}' with data type '{:s}' of "
119 "size {:d} to a different sized type '{:s}' of size {:d}.",
120 array.GetName(), array.GetDataTypeAsString(),
121 array.GetDataTypeSize(), typeid(T).name(), sizeof(T));
122 }
123
124 vtkIdType const nTuples(array.GetNumberOfTuples());
125 int const nComponents(array.GetNumberOfComponents());
126 char const* const array_name(array.GetName());
127
128 auto* const vec = properties.createNewPropertyVector<T>(
129 array_name, type, nComponents);
130 if (!vec)
131 {
132 WARN("Array {:s} could not be converted to PropertyVector.",
133 array_name);
134 return;
135 }
136 vec->reserve(nTuples * nComponents);
137 auto* data_array = static_cast<T*>(array.GetVoidPointer(0));
138 std::copy(&data_array[0],
139 &data_array[nTuples * nComponents],
140 std::back_inserter(*vec));
141 return;
142 }
143};
144
145} // end namespace MeshLib
#define OGS_FATAL(...)
Definition Error.h:26
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40
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....
Definition Properties.h:33
PropertyVector< T > * createNewPropertyVector(std::string_view name, MeshItemType mesh_item_type, std::size_t n_components=1)
Definition Properties.h:17
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.
Definition MeshEnums.h:48
Contains the relevant information when storing a geoscientific raster data.
Definition Raster.h:28