From 8aba0f3bb03b9512e8c48902963fe3741ac1e72d Mon Sep 17 00:00:00 2001 From: marioalexis Date: Sat, 20 Jul 2024 20:03:03 -0300 Subject: [PATCH] Fem: Import/export 1D elements from mesh to vtk - fixes #15541 --- src/Mod/Fem/App/FemVTKTools.cpp | 45 ++++++++++++++++++++++++++++++++- src/Mod/Fem/App/PreCompiled.h | 2 ++ 2 files changed, 46 insertions(+), 1 deletion(-) diff --git a/src/Mod/Fem/App/FemVTKTools.cpp b/src/Mod/Fem/App/FemVTKTools.cpp index 77465d008c..4c952e33e6 100644 --- a/src/Mod/Fem/App/FemVTKTools.cpp +++ b/src/Mod/Fem/App/FemVTKTools.cpp @@ -40,9 +40,11 @@ #include #include #include +#include #include #include #include +#include #include #include #include @@ -166,6 +168,13 @@ void FemVTKTools::importVTKMesh(vtkSmartPointer dataset, FemMesh* me std::vector ids; fillMeshElementIds(cell, ids); switch (cell->GetCellType()) { + // 1D edges + case VTK_LINE: // seg2 + meshds->AddEdgeWithID(ids[0], ids[1], iCell + 1); + break; + case VTK_QUADRATIC_EDGE: // seg3 + meshds->AddEdgeWithID(ids[0], ids[1], ids[2], iCell + 1); + break; // 2D faces case VTK_TRIANGLE: // tria3 meshds->AddFaceWithID(ids[0], ids[1], ids[2], iCell + 1); @@ -282,7 +291,7 @@ void FemVTKTools::importVTKMesh(vtkSmartPointer dataset, FemMesh* me // not handled cases default: { Base::Console().Error( - "Only common 2D and 3D Cells are supported in VTK mesh import\n"); + "Only common 1D, 2D and 3D Cells are supported in VTK mesh import\n"); break; } } @@ -330,6 +339,36 @@ FemMesh* FemVTKTools::readVTKMesh(const char* filename, FemMesh* mesh) return mesh; } +void exportFemMeshEdges(vtkSmartPointer grid, + const SMDS_EdgeIteratorPtr& aEdgeIter) +{ + Base::Console().Log(" Start: VTK mesh builder edges.\n"); + + vtkSmartPointer elemArray = vtkSmartPointer::New(); + std::vector types; + + while (aEdgeIter->more()) { + const SMDS_MeshEdge* aEdge = aEdgeIter->next(); + // edge + if (aEdge->GetEntityType() == SMDSEntity_Edge) { + fillVtkArray(elemArray, types, aEdge); + } + // quadratic edge + else if (aEdge->GetEntityType() == SMDSEntity_Quad_Edge) { + fillVtkArray(elemArray, types, aEdge); + } + else { + throw Base::TypeError("Edge not yet supported by FreeCAD's VTK mesh builder\n"); + } + } + + if (elemArray->GetNumberOfCells() > 0) { + grid->SetCells(types.data(), elemArray); + } + + Base::Console().Log(" End: VTK mesh builder edges.\n"); +} + void exportFemMeshFaces(vtkSmartPointer grid, const SMDS_FaceIteratorPtr& aFaceIter) { @@ -450,6 +489,10 @@ void FemVTKTools::exportVTKMesh(const FemMesh* mesh, Base::Console().Log(" Size of nodes in VTK grid: %i.\n", nNodes); Base::Console().Log(" End: VTK mesh builder nodes.\n"); + // edges + SMDS_EdgeIteratorPtr aEdgeIter = meshDS->edgesIterator(); + exportFemMeshEdges(grid, aEdgeIter); + // faces SMDS_FaceIteratorPtr aFaceIter = meshDS->facesIterator(); exportFemMeshFaces(grid, aFaceIter); diff --git a/src/Mod/Fem/App/PreCompiled.h b/src/Mod/Fem/App/PreCompiled.h index 5d0cc033eb..c075ee3cb8 100644 --- a/src/Mod/Fem/App/PreCompiled.h +++ b/src/Mod/Fem/App/PreCompiled.h @@ -163,12 +163,14 @@ #include #include #include +#include #include #include #include #include #include #include +#include #include #include #include