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/Location.h"
23#include "MeshLib/MeshEnums.h"
24#include "MeshLib/Properties.h"
26
27class vtkUnstructuredGrid; // For conversion vom vtk to ogs mesh
28class vtkDataArray; // For node/cell properties
29
30namespace MeshLib
31{
32class Mesh;
33class Element;
34class Node;
35
40{
41public:
44 vtkUnstructuredGrid* grid,
45 bool const compute_element_neighbors = false,
46 std::string const& mesh_name = "vtkUnstructuredGrid");
47
48private:
49 static void convertScalarArrays(vtkUnstructuredGrid& grid,
50 MeshLib::Mesh& mesh);
51
52 static std::vector<MeshLib::Node*> createNodeVector(
53 std::vector<double> const& elevation,
54 std::vector<int>& node_idx_map,
55 GeoLib::RasterHeader const& header,
56 bool use_elevation);
57
59 static std::vector<MeshLib::Element*> createElementVector(
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,
66 MeshElemType elem_type);
67
69 template <typename T>
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,
75 MeshElemType elem_type)
76 {
77 for (std::size_t i = 0; i < imgHeight; i++)
78 {
79 for (std::size_t j = 0; j < imgWidth; j++)
80 {
81 if (!pix_vis[i * imgWidth + j])
82 {
83 continue;
84 }
85 T val(static_cast<T>(pix_val[i * (imgWidth + 1) + j]));
86 if (elem_type == MeshElemType::TRIANGLE)
87 {
88 prop_vec.push_back(val);
89 prop_vec.push_back(val);
90 }
91 else if (elem_type == MeshElemType::QUAD)
92 {
93 prop_vec.push_back(val);
94 }
95 }
96 }
97 }
98
99 static void convertArray(vtkDataArray& array,
100 MeshLib::Properties& properties,
102
103 template <typename T>
104 static void convertTypedArray(vtkDataArray& array,
105 MeshLib::Properties& properties,
107 {
108 // Keep for debugging purposes, since the vtk array type and the end
109 // type T are not always the same.
110 // DBUG(
111 // "Converting a vtk array '{:s}' with data type '{:s}' of "
112 // "size {:d} to a type '{:s}' of size {:d}.",
113 // array.GetName(), array.GetDataTypeAsString(),
114 // array.GetDataTypeSize(), typeid(T).name(), sizeof(T));
115
116 if (sizeof(T) != array.GetDataTypeSize())
117 {
118 OGS_FATAL(
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));
123 }
124
125 vtkIdType const nTuples(array.GetNumberOfTuples());
126 int const nComponents(array.GetNumberOfComponents());
127 char const* const array_name(array.GetName());
128
129 auto* const vec = properties.createNewPropertyVector<T>(
130 array_name, type, nComponents);
131 if (!vec)
132 {
133 WARN("Array {:s} could not be converted to PropertyVector.",
134 array_name);
135 return;
136 }
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));
142 return;
143 }
144};
145
146} // 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:36
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)
MeshItemType
Definition Location.h:21
MeshElemType
Types of mesh elements supported by OpenGeoSys. Values are from VTKCellType enum.
Definition MeshEnums.h:27
Contains the relevant information when storing a geoscientific raster data.
Definition Raster.h:28