22 #include <vtkBitArray.h>
23 #include <vtkCharArray.h>
24 #include <vtkDoubleArray.h>
25 #include <vtkFloatArray.h>
26 #include <vtkImageData.h>
27 #include <vtkIntArray.h>
28 #include <vtkLongArray.h>
29 #include <vtkLongLongArray.h>
30 #include <vtkPointData.h>
31 #include <vtkSmartPointer.h>
32 #include <vtkUnsignedCharArray.h>
33 #include <vtkUnsignedLongArray.h>
34 #include <vtkUnsignedLongLongArray.h>
38 #include <vtkCellData.h>
39 #include <vtkUnsignedIntArray.h>
40 #include <vtkUnstructuredGrid.h>
46 template <
class T_ELEMENT>
48 const std::vector<MeshLib::Node*>& nodes, vtkIdList*
const node_ids,
49 std::size_t
const element_id)
51 auto** ele_nodes =
new MeshLib::Node*[T_ELEMENT::n_all_nodes];
52 for (
unsigned k(0); k < T_ELEMENT::n_all_nodes; k++)
54 ele_nodes[k] = nodes[node_ids->GetId(k)];
56 return new T_ELEMENT(ele_nodes, element_id);
61 vtkUnstructuredGrid* grid, std::string
const& mesh_name)
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++)
74 coords = grid->GetPoints()->GetPoint(i);
75 nodes[i] =
new MeshLib::Node(coords[0], coords[1], coords[2], i);
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++)
85 grid->GetCellPoints(i, node_ids);
87 int cell_type = grid->GetCellType(i);
92 elem = detail::createElementWithSameNodeOrder<MeshLib::Point>(
98 elem = detail::createElementWithSameNodeOrder<MeshLib::Line>(
104 elem = detail::createElementWithSameNodeOrder<MeshLib::Tri>(
110 elem = detail::createElementWithSameNodeOrder<MeshLib::Quad>(
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)];
126 elem = detail::createElementWithSameNodeOrder<MeshLib::Tet>(
132 elem = detail::createElementWithSameNodeOrder<MeshLib::Hex>(
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)];
152 elem = detail::createElementWithSameNodeOrder<MeshLib::Pyramid>(
159 for (
unsigned j = 0; j < 3; ++j)
161 prism_nodes[j] = nodes[node_ids->GetId(j + 3)];
162 prism_nodes[j + 3] = nodes[node_ids->GetId(j)];
167 case VTK_QUADRATIC_EDGE:
169 elem = detail::createElementWithSameNodeOrder<MeshLib::Line3>(
173 case VTK_QUADRATIC_TRIANGLE:
175 elem = detail::createElementWithSameNodeOrder<MeshLib::Tri6>(
179 case VTK_QUADRATIC_QUAD:
181 elem = detail::createElementWithSameNodeOrder<MeshLib::Quad8>(
185 case VTK_BIQUADRATIC_QUAD:
187 elem = detail::createElementWithSameNodeOrder<MeshLib::Quad9>(
191 case VTK_QUADRATIC_TETRA:
193 elem = detail::createElementWithSameNodeOrder<MeshLib::Tet10>(
197 case VTK_QUADRATIC_HEXAHEDRON:
199 elem = detail::createElementWithSameNodeOrder<MeshLib::Hex20>(
203 case VTK_QUADRATIC_PYRAMID:
206 detail::createElementWithSameNodeOrder<MeshLib::Pyramid13>(
210 case VTK_QUADRATIC_WEDGE:
213 for (
unsigned j = 0; j < 3; ++j)
215 prism_nodes[j] = nodes[node_ids->GetId(j + 3)];
216 prism_nodes[j + 3] = nodes[node_ids->GetId(j)];
218 for (
unsigned j = 0; j < 3; ++j)
220 prism_nodes[6 + j] = nodes[node_ids->GetId(8 - j)];
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)
227 prism_nodes[12 + j] = nodes[node_ids->GetId(11 - j)];
233 ERR(
"VtkMeshConverter::convertUnstructuredGrid(): Unknown mesh "
234 "element type '{:d}'.",
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)
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)
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)
277 *vtkDataArray::SafeDownCast(field_data->GetAbstractArray(i)),
287 if (vtkDoubleArray::SafeDownCast(&array))
289 VtkMeshConverter::convertTypedArray<double>(array, properties, type);
293 if (vtkFloatArray::SafeDownCast(&array))
295 VtkMeshConverter::convertTypedArray<float>(array, properties, type);
299 if (vtkIntArray::SafeDownCast(&array))
301 VtkMeshConverter::convertTypedArray<int>(array, properties, type);
305 if (vtkBitArray::SafeDownCast(&array))
307 VtkMeshConverter::convertTypedArray<bool>(array, properties, type);
311 if (vtkCharArray::SafeDownCast(&array))
313 VtkMeshConverter::convertTypedArray<char>(array, properties, type);
317 if (vtkLongArray::SafeDownCast(&array))
319 VtkMeshConverter::convertTypedArray<long>(array, properties, type);
323 if (vtkLongLongArray::SafeDownCast(&array))
325 VtkMeshConverter::convertTypedArray<long long>(array, properties, type);
329 if (vtkUnsignedLongArray::SafeDownCast(&array))
331 VtkMeshConverter::convertTypedArray<unsigned long>(array, properties,
336 if (vtkUnsignedLongLongArray::SafeDownCast(&array))
338 VtkMeshConverter::convertTypedArray<unsigned long long>(
339 array, properties, type);
343 if (vtkUnsignedIntArray::SafeDownCast(&array))
346 if (std::strncmp(array.GetName(),
"MaterialIDs", 11) == 0)
348 VtkMeshConverter::convertTypedArray<int>(array, properties, type);
352 VtkMeshConverter::convertTypedArray<unsigned>(array, properties,
360 "Array '{:s}' in VTU file uses unsupported data type '{:s}'. The data "
361 "array will not be available.",
362 array.GetName(), array.GetDataTypeAsString());
void ERR(char const *fmt, Args const &... args)
void WARN(char const *fmt, Args const &... args)
Definition of the Mesh class.
Definition of the Node class.
Definition of the VtkMeshConverter class.
Properties & getProperties()
Property manager on mesh items. Class Properties manages scalar, vector or matrix properties....
static MeshLib::Mesh * convertUnstructuredGrid(vtkUnstructuredGrid *grid, std::string const &mesh_name="vtkUnstructuredGrid")
Converts a vtkUnstructuredGrid object to a Mesh.
static void convertScalarArrays(vtkUnstructuredGrid &grid, MeshLib::Mesh &mesh)
static void convertArray(vtkDataArray &array, MeshLib::Properties &properties, MeshLib::MeshItemType type)
MeshLib::Element * createElementWithSameNodeOrder(const std::vector< MeshLib::Node * > &nodes, vtkIdList *const node_ids, std::size_t const element_id)
TemplateElement< MeshLib::PrismRule6 > Prism
TemplateElement< MeshLib::PrismRule15 > Prism15
TemplateElement< MeshLib::HexRule8 > Hex
TemplateElement< MeshLib::QuadRule4 > Quad