// Copyright (c) 2015 GeometryFactory (France). All rights reserved. // // This file is part of CGAL (www.cgal.org). // // $URL$ // $Id$ // SPDX-License-Identifier: GPL-3.0-or-later OR LicenseRef-Commercial // // Author(s) : Sebastien Loriot , // Jane Tournois // #include #include #include "Scene_surface_mesh_item.h" #include "Scene_c3t3_item.h" #include "Scene_polylines_item.h" #include "Scene_points_with_normal_item.h" #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include #include typedef Scene_surface_mesh_item Scene_facegraph_item; typedef Scene_facegraph_item::Face_graph FaceGraph; typedef boost::property_traits::type>::value_type Point; namespace CGAL{ template void extract_segments_from_vtkPointSet(vtkPointSet* poly_data, std::vector< std::vector >& segments) { // get nb of points and cells vtkIdType nb_points = poly_data->GetNumberOfPoints(); vtkIdType nb_cells = poly_data->GetNumberOfCells(); //extract points std::vector point_map(nb_points); for (vtkIdType i = 0; iGetPoint(i, coords); point_map[i]=Point_3(coords[0], coords[1], coords[2]); } //extract segments for (vtkIdType i = 0; iGetCellType(i) != 3 && poly_data->GetCellType(i) != 4) continue; vtkCell* cell_ptr = poly_data->GetCell(i); vtkIdType nb_vertices = cell_ptr->GetNumberOfPoints(); segments.push_back( std::vector() ); for(int j = 0; j < nb_vertices; ++j) { segments.back().push_back(point_map[cell_ptr->GetPointId(j)]); } } } template void polygon_mesh_to_vtkUnstructured(const PM& pmesh,//PolygonMesh const char* filename) { typedef typename boost::graph_traits::vertex_descriptor vertex_descriptor; typedef typename boost::graph_traits::face_descriptor face_descriptor; typedef typename boost::graph_traits::halfedge_descriptor halfedge_descriptor; typedef typename boost::property_map::const_type VPMap; typedef typename boost::property_map_value::type Point_3; VPMap vpmap = get(CGAL::vertex_point, pmesh); vtkPoints* const vtk_points = vtkPoints::New(); vtkCellArray* const vtk_cells = vtkCellArray::New(); vtk_points->Allocate(num_vertices(pmesh)); vtk_cells->Allocate(num_faces(pmesh)); std::map Vids; vtkIdType inum = 0; for(vertex_descriptor v : vertices(pmesh)) { const Point_3& p = get(vpmap, v); vtk_points->InsertNextPoint(CGAL::to_double(p.x()), CGAL::to_double(p.y()), CGAL::to_double(p.z())); Vids[v] = inum++; } for(face_descriptor f : faces(pmesh)) { vtkIdList* cell = vtkIdList::New(); for(halfedge_descriptor h : halfedges_around_face(halfedge(f, pmesh), pmesh)) { cell->InsertNextId(Vids[target(h, pmesh)]); } vtk_cells->InsertNextCell(cell); cell->Delete(); } vtkSmartPointer usg = vtkSmartPointer::New(); usg->SetPoints(vtk_points); vtk_points->Delete(); usg->SetCells(5,vtk_cells); vtk_cells->Delete(); // Combine the two data sets //vtkSmartPointer appendFilter = // vtkSmartPointer::New(); //appendFilter->AddInputData(polydata); //appendFilter->Update(); //vtkSmartPointer unstructuredGrid = // vtkSmartPointer::New(); //unstructuredGrid->ShallowCopy(appendFilter->GetOutput()); // Write the unstructured grid vtkSmartPointer writer = vtkSmartPointer::New(); writer->SetFileName(filename); writer->SetInputData(usg); writer->Write(); } }//end namespace CGAL class Polyhedron_demo_vtk_plugin : public QObject, public CGAL::Three::Polyhedron_demo_io_plugin_interface { Q_OBJECT Q_INTERFACES(CGAL::Three::Polyhedron_demo_io_plugin_interface) Q_PLUGIN_METADATA(IID "com.geometryfactory.PolyhedronDemo.IOPluginInterface/1.90" FILE "vtk_io_plugin.json") public: typedef boost::graph_traits::vertex_descriptor vertex_descriptor; typedef boost::graph_traits::face_descriptor face_descriptor; QString nameFilters() const { return "VTK XML UnstructuredGrid (*.vtu);;VTK PolyData files (*.vtk);; VTK XML PolyData (*.vtp)"; } QString name() const { return "vtk_plugin"; } bool canSave(const CGAL::Three::Scene_item* item) { return (qobject_cast(item) || qobject_cast(item)); } bool save(QFileInfo fileinfo,QList& items) { Scene_item* item = items.front(); std::string extension = fileinfo.suffix().toLower().toStdString(); if ( extension != "vtk" && extension != "vtp" && extension != "vtu") return false; std::string output_filename = fileinfo.absoluteFilePath().toStdString(); const Scene_facegraph_item* poly_item = qobject_cast(item); if (poly_item) { if (extension != "vtp") { if(!CGAL::is_triangle_mesh(*poly_item->polyhedron())) { QMessageBox::warning(0, "Error", "Cannot save a mesh in vtu format if " "it is not pure triangle."); return false; } CGAL::polygon_mesh_to_vtkUnstructured( *poly_item->polyhedron(), output_filename.data()); } else { const FaceGraph* mesh = poly_item->face_graph(); std::ofstream os(output_filename.data()); os << std::setprecision(16); //write header CGAL::IO::write_VTP(os, *mesh); } } else { const Scene_c3t3_item* c3t3_item = qobject_cast(item); if(!c3t3_item || extension != "vtu") return false; std::ofstream os(output_filename.data()); os << std::setprecision(16); const C3t3& c3t3 = c3t3_item->c3t3(); CGAL::IO::output_to_vtu(os, c3t3); } items.pop_front(); return true; } bool canLoad(QFileInfo) const { return true; } QList load(QFileInfo fileinfo, bool& ok, bool add_to_scene) { std::string extension=fileinfo.suffix().toLower().toStdString(); if (extension != "vtk" && extension != "vtp" && extension != "vtu") { ok = false; return QList(); } std::string fname = fileinfo.absoluteFilePath().toStdString(); // Try to read .vtk in a facegraph if(fileinfo.size() == 0) { CGAL::Three::Three::warning( tr("The file you are trying to load is empty.")); Scene_facegraph_item* item = new Scene_facegraph_item(); item->setName(fileinfo.completeBaseName()); ok = true; if(add_to_scene) CGAL::Three::Three::scene()->addItem(item); return QList()< data; vtkSmartPointer obs = vtkSmartPointer::New(); if (extension=="vtp") data = CGAL::IO::internal::read_vtk_file(fname,obs) ->GetOutput(); else if (extension=="vtu") data = CGAL::IO::internal::read_vtk_file(fname,obs) ->GetOutput(); else{ //read non-XML data vtkSmartPointer reader = CGAL::IO::internal::read_vtk_file(fname,obs); data = vtkPolyData::SafeDownCast(reader->GetOutput()); if (!data) data = vtkUnstructuredGrid::SafeDownCast(reader->GetOutput()); } if (obs->GetError()) { QMessageBox msgBox; msgBox.setText("This type of data can't be opened"); msgBox.setInformativeText(QString("VTK error message :\n") .append(QString(obs->GetErrorMessage().data()))); msgBox.setStandardButtons(QMessageBox::Ok); msgBox.setIcon(QMessageBox::Critical); msgBox.exec(); ok = false; return QList(); } if (obs->GetWarning()) { QMessageBox msgBox; msgBox.setText("This file generates a warning"); msgBox.setInformativeText(QString("VTK warning message :\n") .append(QString(obs->GetWarningMessage().data()))); msgBox.setStandardButtons(QMessageBox::Ok); msgBox.setIcon(QMessageBox::Warning); msgBox.exec(); } //check celltypes bool is_polygon_mesh(false), is_c3t3(false), is_polyline(false); for(int i = 0; i< data->GetNumberOfCells(); ++i) { int t = data->GetCellType(i); if( t == 5 || t == 7 || t == 9) //tri, quad or polygon is_polygon_mesh = true; else if(t == 10) //tetrahedron is_c3t3 = true; else if( t == 3 || t == 4) //line or polyline is_polyline = true; } Scene_group_item* group = nullptr; if((is_polygon_mesh && is_c3t3) || (is_polygon_mesh && is_polyline) || (is_c3t3 && is_polyline) ) { group = new Scene_group_item(fileinfo.baseName()); } if(is_polygon_mesh) { auto poly = std::make_unique(); Scene_item* poly_item = nullptr; if (CGAL::IO::internal::vtkPointSet_to_polygon_mesh(data, *poly, CGAL::parameters::default_values())) { poly_item = new Scene_facegraph_item(poly.release()); } else { poly.reset(nullptr); std::vector points; std::vector > polygons; if (CGAL::IO::internal::vtkPointSet_to_polygon_soup(data, points, polygons, CGAL::parameters::default_values())) { auto soup_item = new Scene_polygon_soup_item(); soup_item->load(points, polygons); poly_item = soup_item; } } if(poly_item) { if(group) { poly_item->setName(QString("%1_faces").arg(fileinfo.baseName())); CGAL::Three::Three::scene()->addItem(poly_item); CGAL::Three::Three::scene()->changeGroup(poly_item, group); } else{ poly_item->setName(fileinfo.baseName()); ok = true; if(add_to_scene) CGAL::Three::Three::scene()->addItem(poly_item); return QList()< Facet; // 3 = id typedef std::array Tet; // first 4 = id, fifth = reference Scene_c3t3_item* c3t3_item = new Scene_c3t3_item(); c3t3_item->set_valid(false); //build a triangulation from data: std::vector points; vtkPoints* dataP = data->GetPoints();; for(int i = 0; i< data->GetNumberOfPoints(); ++i) { double *p = dataP->GetPoint(i); points.push_back(Tr::Point(p[0],p[1],p[2])); } std::vector finite_cells; std::vector subdomains; bool has_mesh_domain = data->GetCellData()->HasArray("MeshDomain"); vtkDataArray* domains = data->GetCellData()->GetArray("MeshDomain"); for(int i = 0; i< data->GetNumberOfCells(); ++i) { if(data->GetCellType(i) != 10 ) continue; vtkIdList* pids = data->GetCell(i)->GetPointIds(); Tet cell; for(int j = 0; j<4; ++j) cell[j] = pids->GetId(j); finite_cells.push_back(cell); const auto si = has_mesh_domain ? static_cast(domains->GetComponent(i, 0)) : 1; subdomains.push_back(si); } std::map border_facets; //Preprocessing for build_triangulation //check for orientation and swap in cells if not good. for(std::size_t i=0; ic3t3().triangulation(), points, finite_cells, subdomains, border_facets, false, false, true); for( C3t3::Triangulation::Finite_cells_iterator cit = c3t3_item->c3t3().triangulation().finite_cells_begin(); cit != c3t3_item->c3t3().triangulation().finite_cells_end(); ++cit) { CGAL_assertion(cit->info() >= 0); c3t3_item->c3t3().add_to_complex(cit, cit->info()); for(int i=0; i < 4; ++i) { if(cit->surface_patch_index(i)>0) { c3t3_item->c3t3().add_to_complex(cit, i, cit->surface_patch_index(i)); } } } //if there is no facet in the complex, we add the border facets. if(c3t3_item->c3t3().number_of_facets_in_complex() == 0) { for( C3t3::Triangulation::Finite_facets_iterator fit = c3t3_item->c3t3().triangulation().finite_facets_begin(); fit != c3t3_item->c3t3().triangulation().finite_facets_end(); ++fit) { typedef C3t3::Triangulation::Cell_handle Cell_handle; Cell_handle c = fit->first; Cell_handle nc = c->neighbor(fit->second); // By definition, Subdomain_index() is supposed to be the id of the exterior if(c->subdomain_index() != C3t3::Triangulation::Cell::Subdomain_index() && nc->subdomain_index() == C3t3::Triangulation::Cell::Subdomain_index()) { // Color the border facet with the index of its cell c3t3_item->c3t3().add_to_complex(c, fit->second, c->subdomain_index()); } } } c3t3_item->c3t3_changed(); c3t3_item->resetCutPlane(); if(group) { c3t3_item->setName(QString("%1_tetrahedra").arg(fileinfo.baseName())); CGAL::Three::Three::scene()->addItem(c3t3_item); CGAL::Three::Three::scene()->changeGroup(c3t3_item, group); } else{ c3t3_item->setName(fileinfo.baseName()); ok = true; if(add_to_scene) CGAL::Three::Three::scene()->addItem(c3t3_item); return QList()< > segments; extract_segments_from_vtkPointSet(data,segments); Scene_polylines_item* polyline_item = new Scene_polylines_item(); for(const std::vector& segment : segments) polyline_item->polylines.push_back(segment); if(group) { polyline_item->setName(QString("%1_lines").arg(fileinfo.baseName())); CGAL::Three::Three::scene()->addItem(polyline_item); CGAL::Three::Three::scene()->changeGroup(polyline_item, group); } else{ polyline_item->setName(fileinfo.baseName()); ok = true; if(add_to_scene) CGAL::Three::Three::scene()->addItem(polyline_item); return QList()<addItem(group); return QList()<GetNumberOfPoints(); ++i) { double* p = data->GetPoint(i); point_item->point_set()->insert(Point_3(p[0], p[1], p[2])); } point_item->setName(fileinfo.baseName()); ok = true; if(add_to_scene) CGAL::Three::Three::scene()->addItem(point_item); return QList()<