OGS
MeshLib::VtkMeshConverter Class Reference

Detailed Description

Converter for VtkUnstructured Grids to OGS meshes.

Definition at line 38 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 277 of file VtkMeshConverter.cpp.

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

References convertTypedArray(), OGS_FATAL, and WARN().

Referenced by convertScalarArrays().

◆ convertScalarArrays()

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

Definition at line 242 of file VtkMeshConverter.cpp.

244{
245 vtkPointData* point_data = grid.GetPointData();
246 auto const n_point_arrays =
247 static_cast<unsigned>(point_data->GetNumberOfArrays());
248 for (unsigned i = 0; i < n_point_arrays; ++i)
249 {
250 convertArray(*point_data->GetArray(i),
251 mesh.getProperties(),
253 }
254
255 vtkCellData* cell_data = grid.GetCellData();
256 auto const n_cell_arrays =
257 static_cast<unsigned>(cell_data->GetNumberOfArrays());
258 for (unsigned i = 0; i < n_cell_arrays; ++i)
259 {
260 convertArray(*cell_data->GetArray(i),
261 mesh.getProperties(),
263 }
264
265 vtkFieldData* field_data = grid.GetFieldData();
266 auto const n_field_arrays =
267 static_cast<unsigned>(field_data->GetNumberOfArrays());
268 for (unsigned i = 0; i < n_field_arrays; ++i)
269 {
271 *vtkDataArray::SafeDownCast(field_data->GetAbstractArray(i)),
272 mesh.getProperties(),
274 }
275}
Properties & getProperties()
Definition Mesh.h:136
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 // Keep for debugging purposes, since the vtk array type and the end
108 // type T are not always the same.
109 // DBUG(
110 // "Converting a vtk array '{:s}' with data type '{:s}' of "
111 // "size {:d} to a type '{:s}' of size {:d}.",
112 // array.GetName(), array.GetDataTypeAsString(),
113 // array.GetDataTypeSize(), typeid(T).name(), sizeof(T));
114
115 if (sizeof(T) != array.GetDataTypeSize())
116 {
117 OGS_FATAL(
118 "Trying to convert a vtk array '{:s}' with data type '{:s}' of "
119 "size {:d} to a different sized type '{:s}' of size {:d}.",
120 array.GetName(), array.GetDataTypeAsString(),
121 array.GetDataTypeSize(), typeid(T).name(), sizeof(T));
122 }
123
124 vtkIdType const nTuples(array.GetNumberOfTuples());
125 int const nComponents(array.GetNumberOfComponents());
126 char const* const array_name(array.GetName());
127
128 auto* const vec = properties.createNewPropertyVector<T>(
129 array_name, type, nComponents);
130 if (!vec)
131 {
132 WARN("Array {:s} could not be converted to PropertyVector.",
133 array_name);
134 return;
135 }
136 vec->reserve(nTuples * nComponents);
137 auto* data_array = static_cast<T*>(array.GetVoidPointer(0));
138 std::copy(&data_array[0],
139 &data_array[nTuples * nComponents],
140 std::back_inserter(*vec));
141 return;
142 }
PropertyVector< T > * createNewPropertyVector(std::string_view name, MeshItemType mesh_item_type, std::size_t n_components=1)
Definition Properties.h:17

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

Referenced by convertArray().

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

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

Referenced by ApplicationUtils::computeNaturalCoords(), 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: