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, bool const compute_element_neighbors=false, std::string const &mesh_name="vtkUnstructuredGrid")
 Converts a vtkUnstructuredGrid object to a Mesh.
 

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.
 
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.
 
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 275 of file VtkMeshConverter.cpp.

278{
279 if (vtkDoubleArray::SafeDownCast(&array))
280 {
281 VtkMeshConverter::convertTypedArray<double>(array, properties, type);
282 return;
283 }
284
285 if (vtkFloatArray::SafeDownCast(&array))
286 {
287 VtkMeshConverter::convertTypedArray<float>(array, properties, type);
288 return;
289 }
290
291 if (vtkBitArray::SafeDownCast(&array))
292 {
293 VtkMeshConverter::convertTypedArray<bool>(array, properties, type);
294 return;
295 }
296
297 // This also covers the vtkTypeInt8Array type, which is derived from the
298 // vtkCharArray type.
299 if (vtkCharArray::SafeDownCast(&array) ||
300 vtkSignedCharArray::SafeDownCast(&array))
301 {
302 if constexpr (sizeof(std::int8_t) != sizeof(char))
303 {
304 OGS_FATAL(
305 "Array '{:s}' in VTU file uses unsupported data type '{:s}' "
306 "not convertible to char.",
307 array.GetName(), array.GetDataTypeAsString());
308 }
309 VtkMeshConverter::convertTypedArray<char>(array, properties, type);
310 return;
311 }
312
313 // This also covers the vtkTypeInt16Array type, which is derived from the
314 // vtkShortArray type.
315 if (vtkShortArray::SafeDownCast(&array))
316 {
317 if constexpr (sizeof(std::int16_t) != sizeof(short))
318 {
319 OGS_FATAL(
320 "Array '{:s}' in VTU file uses unsupported data type '{:s}' "
321 "not convertible to short.",
322 array.GetName(), array.GetDataTypeAsString());
323 }
324 VtkMeshConverter::convertTypedArray<short>(array, properties, type);
325 return;
326 }
327
328 // This also covers the vtkTypeInt32Array type, which is derived from the
329 // vtkIntArray type.
330 if (vtkIntArray::SafeDownCast(&array))
331 {
332 if constexpr (sizeof(std::int32_t) != sizeof(int))
333 {
334 OGS_FATAL(
335 "Array '{:s}' in VTU file uses unsupported data type '{:s}' "
336 "not convertible to int.",
337 array.GetName(), array.GetDataTypeAsString());
338 }
339 VtkMeshConverter::convertTypedArray<int>(array, properties, type);
340 return;
341 }
342
343 // This is not a unique conversion depending on the sizes of int, long, and
344 // long long. Converting to the apparently smaller type.
345 if (vtkLongArray::SafeDownCast(&array))
346 {
347 if constexpr (sizeof(long) == sizeof(int))
348 {
349 VtkMeshConverter::convertTypedArray<int>(array, properties, type);
350 return;
351 }
352 if constexpr (sizeof(long long) == sizeof(long))
353 {
354 VtkMeshConverter::convertTypedArray<long>(array, properties, type);
355 return;
356 }
357
358 OGS_FATAL(
359 "Array '{:s}' in VTU file uses unsupported data type '{:s}' ({:d}) "
360 "not convertible to int ({:d}), or long ({:d}) types.",
361 array.GetName(), array.GetDataTypeAsString(),
362 array.GetDataTypeSize(), sizeof(int), sizeof(long));
363 }
364
365 // Ensure, vtkTypeInt64Array and vtkLongLongArray are using the same type.
366 static_assert(sizeof(std::int64_t) == sizeof(long long));
367
368 // This also covers the vtkTypeInt64Array type, which is derived from the
369 // vtkLongLongArray type.
370 // Converting to the apparently smaller type.
371 if (vtkLongLongArray::SafeDownCast(&array))
372 {
373 if constexpr (sizeof(long long) == sizeof(long))
374 {
375 VtkMeshConverter::convertTypedArray<long>(array, properties, type);
376 return;
377 }
378 else // ll > l. Other cases are not possible because ll >= l in c++.
379 {
380 VtkMeshConverter::convertTypedArray<long long>(array, properties,
381 type);
382 return;
383 }
384 }
385
386 // This also covers the vtkTypeUInt8Array type, which is derived from the
387 // vtkUnsignedCharArray type.
388 if (vtkUnsignedCharArray::SafeDownCast(&array))
389 {
390 if constexpr (sizeof(std::uint8_t) != sizeof(unsigned char))
391 {
392 OGS_FATAL(
393 "Array '{:s}' in VTU file uses unsupported data type '{:s}' "
394 "not convertible to unsigned char.",
395 array.GetName(), array.GetDataTypeAsString());
396 }
397 VtkMeshConverter::convertTypedArray<unsigned char>(array, properties,
398 type);
399 return;
400 }
401
402 // This also covers the vtkTypeUInt16Array type, which is derived from the
403 // vtkUnsignedShortArray type.
404 if (vtkUnsignedShortArray::SafeDownCast(&array))
405 {
406 if constexpr (sizeof(std::uint16_t) != sizeof(unsigned short))
407 {
408 OGS_FATAL(
409 "Array '{:s}' in VTU file uses unsupported data type '{:s}' "
410 "not convertible to unsigned short.",
411 array.GetName(), array.GetDataTypeAsString());
412 }
413 VtkMeshConverter::convertTypedArray<unsigned short>(array, properties,
414 type);
415 return;
416 }
417
418 // This also covers the vtkTypeUInt32Array type, which is derived from the
419 // vtkUnsignedIntArray type.
420 if (vtkUnsignedIntArray::SafeDownCast(&array))
421 {
422 if constexpr (sizeof(std::uint32_t) != sizeof(unsigned))
423 {
424 OGS_FATAL(
425 "Array '{:s}' in VTU file uses unsupported data type '{:s}' "
426 "not convertible to unsigned.",
427 array.GetName(), array.GetDataTypeAsString());
428 }
429
430 // MaterialIDs are assumed to be integers
431 if (std::strncmp(array.GetName(), "MaterialIDs", 11) == 0)
432 {
433 VtkMeshConverter::convertTypedArray<int>(array, properties, type);
434 }
435 else
436 {
437 VtkMeshConverter::convertTypedArray<unsigned>(array, properties,
438 type);
439 }
440
441 return;
442 }
443
444 // This is not a unique conversion depending on the sizes of unsigned,
445 // unsigned long, and unsigned long long. Converting to the apparently
446 // smaller type.
447 if (vtkUnsignedLongArray::SafeDownCast(&array))
448 {
449 if constexpr (sizeof(unsigned long) == sizeof(unsigned))
450 {
451 VtkMeshConverter::convertTypedArray<unsigned>(array, properties,
452 type);
453 return;
454 }
455 if constexpr (sizeof(unsigned long long) == sizeof(unsigned long))
456 {
457 VtkMeshConverter::convertTypedArray<unsigned long>(
458 array, properties, type);
459 return;
460 }
461
462 OGS_FATAL(
463 "Array '{:s}' in VTU file uses unsupported data type '{:s}' ({:d}) "
464 "not convertible to unsigned ({:d}), or unsigned long ({:d}) "
465 "types.",
466 array.GetName(), array.GetDataTypeAsString(),
467 array.GetDataTypeSize(), sizeof(unsigned), sizeof(unsigned long));
468 }
469
470 // Ensure, vtkTypeUInt64Array and vtkUnsignedLongLongArray are using the
471 // same type.
472 static_assert(sizeof(std::uint64_t) == sizeof(unsigned long long));
473
474 // This also covers the vtkTypeUInt64Array type, which is derived from the
475 // vtkUnsignedLongLongArray type.
476 // Converting to the apparently smaller type.
477 if (vtkUnsignedLongLongArray::SafeDownCast(&array))
478 {
479 if constexpr (sizeof(unsigned long long) == sizeof(unsigned long))
480 {
481 VtkMeshConverter::convertTypedArray<unsigned long>(
482 array, properties, type);
483 return;
484 }
485 else // ull > ul. Other cases are not possible because ull >= ul in
486 // c++.
487 {
488 VtkMeshConverter::convertTypedArray<unsigned long long>(
489 array, properties, type);
490 return;
491 }
492 }
493
494 WARN(
495 "Array '{:s}' in VTU file uses unsupported data type '{:s}' of size "
496 "{:d}. The data array will not be available.",
497 array.GetName(), array.GetDataTypeAsString(), array.GetDataTypeSize());
498}
#define OGS_FATAL(...)
Definition Error.h:26
void WARN(fmt::format_string< Args... > fmt, Args &&... args)
Definition Logging.h:40

References OGS_FATAL, and WARN().

Referenced by convertScalarArrays().

◆ convertScalarArrays()

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

Definition at line 240 of file VtkMeshConverter.cpp.

242{
243 vtkPointData* point_data = grid.GetPointData();
244 auto const n_point_arrays =
245 static_cast<unsigned>(point_data->GetNumberOfArrays());
246 for (unsigned i = 0; i < n_point_arrays; ++i)
247 {
248 convertArray(*point_data->GetArray(i),
249 mesh.getProperties(),
251 }
252
253 vtkCellData* cell_data = grid.GetCellData();
254 auto const n_cell_arrays =
255 static_cast<unsigned>(cell_data->GetNumberOfArrays());
256 for (unsigned i = 0; i < n_cell_arrays; ++i)
257 {
258 convertArray(*cell_data->GetArray(i),
259 mesh.getProperties(),
261 }
262
263 vtkFieldData* field_data = grid.GetFieldData();
264 auto const n_field_arrays =
265 static_cast<unsigned>(field_data->GetNumberOfArrays());
266 for (unsigned i = 0; i < n_field_arrays; ++i)
267 {
269 *vtkDataArray::SafeDownCast(field_data->GetAbstractArray(i)),
270 mesh.getProperties(),
272 }
273}
Properties & getProperties()
Definition Mesh.h:134
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 104 of file VtkMeshConverter.h.

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 }
PropertyVector< T > * createNewPropertyVector(std::string_view name, MeshItemType mesh_item_type, std::size_t n_components=1)

References MeshLib::Properties::createNewPropertyVector(), OGS_FATAL, and WARN().

◆ convertUnstructuredGrid()

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

Converts a vtkUnstructuredGrid object to a Mesh.

Definition at line 65 of file VtkMeshConverter.cpp.

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

References convertScalarArrays(), ERR(), and MeshLib::Node.

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 70 of file VtkMeshConverter.h.

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 }

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


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