OGS
MeshLib::VtkMeshConverter Class Reference

Detailed Description

Converter for VtkUnstructured Grids to OGS meshes.

Definition at line 39 of file VtkMeshConverter.h.

#include <VtkMeshConverter.h>

Static Public Member Functions

static MeshLib::MeshconvertUnstructuredGrid (vtkUnstructuredGrid *grid, std::string const &mesh_name="vtkUnstructuredGrid")
 Converts a vtkUnstructuredGrid object to a Mesh. More...
 

Static Private Member Functions

static void convertScalarArrays (vtkUnstructuredGrid &grid, MeshLib::Mesh &mesh)
 
static std::vector< MeshLib::Node * > createNodeVector (std::vector< double > const &elevation, std::vector< int > &node_idx_map, GeoLib::RasterHeader const &header, bool use_elevation)
 
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. More...
 
template<typename T >
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. More...
 
static void convertArray (vtkDataArray &array, MeshLib::Properties &properties, MeshLib::MeshItemType type)
 
template<typename T >
static void convertTypedArray (vtkDataArray &array, MeshLib::Properties &properties, MeshLib::MeshItemType type)
 

Member Function Documentation

◆ convertArray()

void MeshLib::VtkMeshConverter::convertArray ( vtkDataArray &  array,
MeshLib::Properties properties,
MeshLib::MeshItemType  type 
)
staticprivate

Definition at line 283 of file VtkMeshConverter.cpp.

286 {
287  if (vtkDoubleArray::SafeDownCast(&array))
288  {
289  VtkMeshConverter::convertTypedArray<double>(array, properties, type);
290  return;
291  }
292 
293  if (vtkFloatArray::SafeDownCast(&array))
294  {
295  VtkMeshConverter::convertTypedArray<float>(array, properties, type);
296  return;
297  }
298 
299  if (vtkIntArray::SafeDownCast(&array))
300  {
301  VtkMeshConverter::convertTypedArray<int>(array, properties, type);
302  return;
303  }
304 
305  if (vtkBitArray::SafeDownCast(&array))
306  {
307  VtkMeshConverter::convertTypedArray<bool>(array, properties, type);
308  return;
309  }
310 
311  if (vtkCharArray::SafeDownCast(&array))
312  {
313  VtkMeshConverter::convertTypedArray<char>(array, properties, type);
314  return;
315  }
316 
317  if (vtkLongArray::SafeDownCast(&array))
318  {
319  VtkMeshConverter::convertTypedArray<long>(array, properties, type);
320  return;
321  }
322 
323  if (vtkLongLongArray::SafeDownCast(&array))
324  {
325  VtkMeshConverter::convertTypedArray<long long>(array, properties, type);
326  return;
327  }
328 
329  if (vtkUnsignedLongArray::SafeDownCast(&array))
330  {
331  VtkMeshConverter::convertTypedArray<unsigned long>(array, properties,
332  type);
333  return;
334  }
335 
336  if (vtkUnsignedLongLongArray::SafeDownCast(&array))
337  {
338  VtkMeshConverter::convertTypedArray<unsigned long long>(
339  array, properties, type);
340  return;
341  }
342 
343  if (vtkUnsignedIntArray::SafeDownCast(&array))
344  {
345  // MaterialIDs are assumed to be integers
346  if (std::strncmp(array.GetName(), "MaterialIDs", 11) == 0)
347  {
348  VtkMeshConverter::convertTypedArray<int>(array, properties, type);
349  }
350  else
351  {
352  VtkMeshConverter::convertTypedArray<unsigned>(array, properties,
353  type);
354  }
355 
356  return;
357  }
358 
359  WARN(
360  "Array '{:s}' in VTU file uses unsupported data type '{:s}'. The data "
361  "array will not be available.",
362  array.GetName(), array.GetDataTypeAsString());
363 }
void WARN(char const *fmt, Args const &... args)
Definition: Logging.h:37

References WARN().

Referenced by convertScalarArrays().

◆ convertScalarArrays()

void MeshLib::VtkMeshConverter::convertScalarArrays ( vtkUnstructuredGrid &  grid,
MeshLib::Mesh mesh 
)
staticprivate

Definition at line 248 of file VtkMeshConverter.cpp.

250 {
251  vtkPointData* point_data = grid.GetPointData();
252  auto const n_point_arrays =
253  static_cast<unsigned>(point_data->GetNumberOfArrays());
254  for (unsigned i = 0; i < n_point_arrays; ++i)
255  {
256  convertArray(*point_data->GetArray(i),
257  mesh.getProperties(),
259  }
260 
261  vtkCellData* cell_data = grid.GetCellData();
262  auto const n_cell_arrays =
263  static_cast<unsigned>(cell_data->GetNumberOfArrays());
264  for (unsigned i = 0; i < n_cell_arrays; ++i)
265  {
266  convertArray(*cell_data->GetArray(i),
267  mesh.getProperties(),
269  }
270 
271  vtkFieldData* field_data = grid.GetFieldData();
272  auto const n_field_arrays =
273  static_cast<unsigned>(field_data->GetNumberOfArrays());
274  for (unsigned i = 0; i < n_field_arrays; ++i)
275  {
276  convertArray(
277  *vtkDataArray::SafeDownCast(field_data->GetAbstractArray(i)),
278  mesh.getProperties(),
280  }
281 }
Properties & getProperties()
Definition: Mesh.h:123
static void convertArray(vtkDataArray &array, MeshLib::Properties &properties, MeshLib::MeshItemType type)

References MeshLib::Cell, convertArray(), MeshLib::Mesh::getProperties(), MeshLib::IntegrationPoint, and MeshLib::Node.

Referenced by convertUnstructuredGrid().

◆ convertTypedArray()

template<typename T >
static void MeshLib::VtkMeshConverter::convertTypedArray ( vtkDataArray &  array,
MeshLib::Properties properties,
MeshLib::MeshItemType  type 
)
inlinestaticprivate

Definition at line 103 of file VtkMeshConverter.h.

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  }
PropertyVector< T > * createNewPropertyVector(std::string const &name, MeshItemType mesh_item_type, std::size_t n_components=1)
void copy(PETScVector const &x, PETScVector &y)
Definition: LinAlg.cpp:37

References MathLib::LinAlg::copy(), MeshLib::Properties::createNewPropertyVector(), and WARN().

◆ convertUnstructuredGrid()

MeshLib::Mesh * MeshLib::VtkMeshConverter::convertUnstructuredGrid ( vtkUnstructuredGrid *  grid,
std::string const &  mesh_name = "vtkUnstructuredGrid" 
)
static

Converts a vtkUnstructuredGrid object to a Mesh.

Definition at line 60 of file VtkMeshConverter.cpp.

62 {
63  if (!grid)
64  {
65  return nullptr;
66  }
67 
68  // set mesh nodes
69  const std::size_t nNodes = grid->GetPoints()->GetNumberOfPoints();
70  std::vector<MeshLib::Node*> nodes(nNodes);
71  double* coords = nullptr;
72  for (std::size_t i = 0; i < nNodes; i++)
73  {
74  coords = grid->GetPoints()->GetPoint(i);
75  nodes[i] = new MeshLib::Node(coords[0], coords[1], coords[2], i);
76  }
77 
78  // set mesh elements
79  const std::size_t nElems = grid->GetNumberOfCells();
80  std::vector<MeshLib::Element*> elements(nElems);
81  auto node_ids = vtkSmartPointer<vtkIdList>::New();
82  for (std::size_t i = 0; i < nElems; i++)
83  {
84  MeshLib::Element* elem;
85  grid->GetCellPoints(i, node_ids);
86 
87  int cell_type = grid->GetCellType(i);
88  switch (cell_type)
89  {
90  case VTK_VERTEX:
91  {
92  elem = detail::createElementWithSameNodeOrder<MeshLib::Point>(
93  nodes, node_ids, i);
94  break;
95  }
96  case VTK_LINE:
97  {
98  elem = detail::createElementWithSameNodeOrder<MeshLib::Line>(
99  nodes, node_ids, i);
100  break;
101  }
102  case VTK_TRIANGLE:
103  {
104  elem = detail::createElementWithSameNodeOrder<MeshLib::Tri>(
105  nodes, node_ids, i);
106  break;
107  }
108  case VTK_QUAD:
109  {
110  elem = detail::createElementWithSameNodeOrder<MeshLib::Quad>(
111  nodes, node_ids, i);
112  break;
113  }
114  case VTK_PIXEL:
115  {
116  auto** quad_nodes = new MeshLib::Node*[4];
117  quad_nodes[0] = nodes[node_ids->GetId(0)];
118  quad_nodes[1] = nodes[node_ids->GetId(1)];
119  quad_nodes[2] = nodes[node_ids->GetId(3)];
120  quad_nodes[3] = nodes[node_ids->GetId(2)];
121  elem = new MeshLib::Quad(quad_nodes, i);
122  break;
123  }
124  case VTK_TETRA:
125  {
126  elem = detail::createElementWithSameNodeOrder<MeshLib::Tet>(
127  nodes, node_ids, i);
128  break;
129  }
130  case VTK_HEXAHEDRON:
131  {
132  elem = detail::createElementWithSameNodeOrder<MeshLib::Hex>(
133  nodes, node_ids, i);
134  break;
135  }
136  case VTK_VOXEL:
137  {
138  auto** voxel_nodes = new MeshLib::Node*[8];
139  voxel_nodes[0] = nodes[node_ids->GetId(0)];
140  voxel_nodes[1] = nodes[node_ids->GetId(1)];
141  voxel_nodes[2] = nodes[node_ids->GetId(3)];
142  voxel_nodes[3] = nodes[node_ids->GetId(2)];
143  voxel_nodes[4] = nodes[node_ids->GetId(4)];
144  voxel_nodes[5] = nodes[node_ids->GetId(5)];
145  voxel_nodes[6] = nodes[node_ids->GetId(7)];
146  voxel_nodes[7] = nodes[node_ids->GetId(6)];
147  elem = new MeshLib::Hex(voxel_nodes, i);
148  break;
149  }
150  case VTK_PYRAMID:
151  {
152  elem = detail::createElementWithSameNodeOrder<MeshLib::Pyramid>(
153  nodes, node_ids, i);
154  break;
155  }
156  case VTK_WEDGE:
157  {
158  auto** prism_nodes = new MeshLib::Node*[6];
159  for (unsigned j = 0; j < 3; ++j)
160  {
161  prism_nodes[j] = nodes[node_ids->GetId(j + 3)];
162  prism_nodes[j + 3] = nodes[node_ids->GetId(j)];
163  }
164  elem = new MeshLib::Prism(prism_nodes, i);
165  break;
166  }
167  case VTK_QUADRATIC_EDGE:
168  {
169  elem = detail::createElementWithSameNodeOrder<MeshLib::Line3>(
170  nodes, node_ids, i);
171  break;
172  }
173  case VTK_QUADRATIC_TRIANGLE:
174  {
175  elem = detail::createElementWithSameNodeOrder<MeshLib::Tri6>(
176  nodes, node_ids, i);
177  break;
178  }
179  case VTK_QUADRATIC_QUAD:
180  {
181  elem = detail::createElementWithSameNodeOrder<MeshLib::Quad8>(
182  nodes, node_ids, i);
183  break;
184  }
185  case VTK_BIQUADRATIC_QUAD:
186  {
187  elem = detail::createElementWithSameNodeOrder<MeshLib::Quad9>(
188  nodes, node_ids, i);
189  break;
190  }
191  case VTK_QUADRATIC_TETRA:
192  {
193  elem = detail::createElementWithSameNodeOrder<MeshLib::Tet10>(
194  nodes, node_ids, i);
195  break;
196  }
197  case VTK_QUADRATIC_HEXAHEDRON:
198  {
199  elem = detail::createElementWithSameNodeOrder<MeshLib::Hex20>(
200  nodes, node_ids, i);
201  break;
202  }
203  case VTK_QUADRATIC_PYRAMID:
204  {
205  elem =
206  detail::createElementWithSameNodeOrder<MeshLib::Pyramid13>(
207  nodes, node_ids, i);
208  break;
209  }
210  case VTK_QUADRATIC_WEDGE:
211  {
212  auto** prism_nodes = new MeshLib::Node*[15];
213  for (unsigned j = 0; j < 3; ++j)
214  {
215  prism_nodes[j] = nodes[node_ids->GetId(j + 3)];
216  prism_nodes[j + 3] = nodes[node_ids->GetId(j)];
217  }
218  for (unsigned j = 0; j < 3; ++j)
219  {
220  prism_nodes[6 + j] = nodes[node_ids->GetId(8 - j)];
221  }
222  prism_nodes[9] = nodes[node_ids->GetId(12)];
223  prism_nodes[10] = nodes[node_ids->GetId(14)];
224  prism_nodes[11] = nodes[node_ids->GetId(13)];
225  for (unsigned j = 0; j < 3; ++j)
226  {
227  prism_nodes[12 + j] = nodes[node_ids->GetId(11 - j)];
228  }
229  elem = new MeshLib::Prism15(prism_nodes, i);
230  break;
231  }
232  default:
233  ERR("VtkMeshConverter::convertUnstructuredGrid(): Unknown mesh "
234  "element type '{:d}'.",
235  cell_type);
236  return nullptr;
237  }
238 
239  elements[i] = elem;
240  }
241 
242  MeshLib::Mesh* mesh = new MeshLib::Mesh(mesh_name, nodes, elements);
243  convertScalarArrays(*grid, *mesh);
244 
245  return mesh;
246 }
void ERR(char const *fmt, Args const &... args)
Definition: Logging.h:42
static void convertScalarArrays(vtkUnstructuredGrid &grid, MeshLib::Mesh &mesh)
TemplateElement< MeshLib::PrismRule6 > Prism
Definition: Prism.h:25
TemplateElement< MeshLib::PrismRule15 > Prism15
Definition: Prism.h:26
TemplateElement< MeshLib::HexRule8 > Hex
Definition: Hex.h:25
TemplateElement< MeshLib::QuadRule4 > Quad
Definition: Quad.h:28

References convertScalarArrays(), and ERR().

Referenced by VtkVisPipelineView::convertVTKToOGSMesh(), MeshLib::IO::VtuInterface::readVTKFile(), and MeshLib::IO::VtuInterface::readVTUFile().

◆ createElementVector()

static std::vector<MeshLib::Element*> MeshLib::VtkMeshConverter::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 
)
staticprivate

Creates a mesh element vector based on image data.

◆ createNodeVector()

static std::vector<MeshLib::Node*> MeshLib::VtkMeshConverter::createNodeVector ( std::vector< double > const &  elevation,
std::vector< int > &  node_idx_map,
GeoLib::RasterHeader const &  header,
bool  use_elevation 
)
staticprivate

◆ fillPropertyVector()

template<typename T >
static void MeshLib::VtkMeshConverter::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 
)
inlinestaticprivate

Creates a scalar array/mesh property based on pixel values.

Definition at line 69 of file VtkMeshConverter.h.

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  }

References MeshLib::QUAD, and MeshLib::TRIANGLE.


The documentation for this class was generated from the following files: