OGS
VtkMeshConverter.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#pragma once
5
6#include <vtkDataArray.h>
7#include <vtkType.h>
8
10#include "BaseLib/Logging.h"
11#include "GeoLib/Raster.h"
12#include "MeshLib/MeshEnums.h"
13#include "MeshLib/Properties.h"
15
16class vtkUnstructuredGrid; // For conversion vom vtk to ogs mesh
17class vtkDataArray; // For node/cell properties
18
19namespace MeshLib
20{
21class Mesh;
22class Element;
23class Node;
24
29{
30public:
33 vtkUnstructuredGrid* grid,
34 bool const compute_element_neighbors = false,
35 std::string const& mesh_name = "vtkUnstructuredGrid");
36
37private:
38 static void convertScalarArrays(vtkUnstructuredGrid& grid,
39 MeshLib::Mesh& mesh);
40
41 static std::vector<MeshLib::Node*> createNodeVector(
42 std::vector<double> const& elevation,
43 std::vector<int>& node_idx_map,
44 GeoLib::RasterHeader const& header,
45 bool use_elevation);
46
48 static std::vector<MeshLib::Element*> createElementVector(
49 std::vector<double> const& pix_val,
50 std::vector<bool> const& pix_vis,
51 std::vector<MeshLib::Node*> const& nodes,
52 std::vector<int> const& node_idx_map,
53 std::size_t const imgHeight,
54 std::size_t const imgWidth,
55 MeshElemType elem_type);
56
58 template <typename T>
60 std::vector<double> const& pix_val,
61 std::vector<bool> const& pix_vis,
62 const std::size_t& imgHeight,
63 const std::size_t& imgWidth,
64 MeshElemType elem_type)
65 {
66 for (std::size_t i = 0; i < imgHeight; i++)
67 {
68 for (std::size_t j = 0; j < imgWidth; j++)
69 {
70 if (!pix_vis[i * imgWidth + j])
71 {
72 continue;
73 }
74 T val(static_cast<T>(pix_val[i * (imgWidth + 1) + j]));
75 if (elem_type == MeshElemType::TRIANGLE)
76 {
77 prop_vec.push_back(val);
78 prop_vec.push_back(val);
79 }
80 else if (elem_type == MeshElemType::QUAD)
81 {
82 prop_vec.push_back(val);
83 }
84 }
85 }
86 }
87
88 static void convertArray(vtkDataArray& array,
89 MeshLib::Properties& properties,
91
92 template <typename T>
93 static void convertTypedArray(vtkDataArray& array,
94 MeshLib::Properties& properties,
96 {
97 // Keep for debugging purposes, since the vtk array type and the end
98 // type T are not always the same.
99 // DBUG(
100 // "Converting a vtk array '{:s}' with data type '{:s}' of "
101 // "size {:d} to a type '{:s}' of size {:d}.",
102 // array.GetName(), array.GetDataTypeAsString(),
103 // array.GetDataTypeSize(), typeid(T).name(), sizeof(T));
104
105 if (sizeof(T) != array.GetDataTypeSize())
106 {
107 OGS_FATAL(
108 "Trying to convert a vtk array '{:s}' with data type '{:s}' of "
109 "size {:d} to a different sized type '{:s}' of size {:d}.",
110 array.GetName(), array.GetDataTypeAsString(),
111 array.GetDataTypeSize(), BaseLib::typeToString<T>(), sizeof(T));
112 }
113
114 vtkIdType const nTuples(array.GetNumberOfTuples());
115 int const nComponents(array.GetNumberOfComponents());
116 char const* const array_name(array.GetName());
117
118 auto* const vec = properties.createNewPropertyVector<T>(
119 array_name, type, nTuples, nComponents);
120 if (!vec)
121 {
122 WARN("Array {:s} could not be converted to PropertyVector.",
123 array_name);
124 return;
125 }
126 auto* data_array = static_cast<T*>(array.GetVoidPointer(0));
127 vec->assign(std::span<T>(data_array, nTuples * nComponents));
128 return;
129 }
130};
131
132} // end namespace MeshLib
#define OGS_FATAL(...)
Definition Error.h:19
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:34
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)
constexpr void assign(R &&r)
constexpr void push_back(const PROP_VAL_TYPE &value)
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)
std::string typeToString()
MeshElemType
Types of mesh elements supported by OpenGeoSys. Values are from VTKCellType enum.
Definition MeshEnums.h:37
Contains the relevant information when storing a geoscientific raster data.
Definition Raster.h:18