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"
25 #include "MeshLib/PropertyVector.h"
26 
27 class vtkUnstructuredGrid; // For conversion vom vtk to ogs mesh
28 class vtkDataArray; // For node/cell properties
29 
30 namespace MeshLib
31 {
32 class Mesh;
33 class Element;
34 class Node;
35 
40 {
41 public:
44  vtkUnstructuredGrid* grid,
45  std::string const& mesh_name = "vtkUnstructuredGrid");
46 
47 private:
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,
100  MeshLib::MeshItemType type);
101 
102  template <typename T>
103  static void convertTypedArray(vtkDataArray& array,
104  MeshLib::Properties& properties,
106  {
107  vtkIdType const nTuples(array.GetNumberOfTuples());
108  int const nComponents(array.GetNumberOfComponents());
109  char const* const array_name(array.GetName());
110 
111  auto* const vec = properties.createNewPropertyVector<T>(
112  array_name, type, nComponents);
113  if (!vec)
114  {
115  WARN("Array {:s} could not be converted to PropertyVector.",
116  array_name);
117  return;
118  }
119  vec->reserve(nTuples * nComponents);
120  auto* data_array = static_cast<T*>(array.GetVoidPointer(0));
121  std::copy(&data_array[0],
122  &data_array[nTuples * nComponents],
123  std::back_inserter(*vec));
124  return;
125  }
126 };
127 
128 } // end namespace MeshLib
void WARN(char const *fmt, Args const &... args)
Definition: Logging.h:37
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 const &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 MeshLib::Mesh * convertUnstructuredGrid(vtkUnstructuredGrid *grid, std::string const &mesh_name="vtkUnstructuredGrid")
Converts a vtkUnstructuredGrid object to a Mesh.
static void convertTypedArray(vtkDataArray &array, MeshLib::Properties &properties, MeshLib::MeshItemType type)
static void convertScalarArrays(vtkUnstructuredGrid &grid, MeshLib::Mesh &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)
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37
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:25