Skip to content
New issue

Have a question about this project? Sign up for a free GitHub account to open an issue and contact its maintainers and the community.

By clicking “Sign up for GitHub”, you agree to our terms of service and privacy statement. We’ll occasionally send you account related emails.

Already on GitHub? Sign in to your account

Handle arbitrary face order in polyhedral vtk face mesh #2009

Merged
Merged
Show file tree
Hide file tree
Changes from all commits
Commits
File filter

Filter by extension

Filter by extension

Conversations
Failed to load comments.
Loading
Jump to
Jump to file
Failed to load files.
Loading
Diff view
Diff view
78 changes: 70 additions & 8 deletions arcane/src/arcane/std/VtkPolyhedralMeshIOService.cc
Original file line number Diff line number Diff line change
@@ -1,11 +1,11 @@
// -*- tab-width: 2; indent-tabs-mode: nil; coding: utf-8-with-signature -*-
//-----------------------------------------------------------------------------
// Copyright 2000-2024 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
// Copyright 2000-2025 CEA (www.cea.fr) IFPEN (www.ifpenergiesnouvelles.com)
// See the top-level COPYRIGHT file for details.
// SPDX-License-Identifier: Apache-2.0
//-----------------------------------------------------------------------------
/*---------------------------------------------------------------------------*/
/* VtkPolyhedralMeshIOService (C) 2000-2024 */
/* VtkPolyhedralMeshIOService (C) 2000-2025 */
/* */
/* Read/write fools for polyhedral mesh with vtk file format */
/*---------------------------------------------------------------------------*/
Expand All @@ -17,6 +17,7 @@
#include <numeric>
#include <functional>
#include <memory>
#include <iterator>

#include <vtkUnstructuredGrid.h>
#include <vtkUnstructuredGridReader.h>
Expand All @@ -30,6 +31,8 @@
#include <vtkDataArrayAccessor.h>
#include <vtkPolyDataReader.h>
#include <vtkPolyData.h>
#include <vtkSmartPointer.h>
#include <vtkCellArray.h>

#include <arccore/base/Ref.h>
#include <arccore/base/String.h>
Expand Down Expand Up @@ -123,6 +126,9 @@ class VtkPolyhedralMeshIOService
Int64ConstArrayView faceNodes();
Int32ConstArrayView faceNbNodes();

Int64ConstArrayView faceNodesInFaceMesh();
Int32ConstArrayView faceNbNodesInFaceMesh();

Int64ConstArrayView edgeNodes();
Int32ConstArrayView edgeNbNodes();

Expand Down Expand Up @@ -163,7 +169,7 @@ class VtkPolyhedralMeshIOService
vtkCellData* faceData();

bool isEmpty() const noexcept { return m_is_empty; }
bool doRead() const noexcept { return m_do_read; }
bool doRead() const noexcept { return m_do_read; }

private:

Expand All @@ -188,6 +194,8 @@ class VtkPolyhedralMeshIOService
Int32UniqueArray m_node_nb_cells, m_node_nb_faces, m_node_nb_edges;
Int32UniqueArray m_cell_nb_edges, m_face_nb_edges, m_face_uid_indexes;
Int32UniqueArray m_cell_face_indexes, m_edge_nb_nodes;
Int64UniqueArray m_face_node_uids_in_face_mesh;
Int32UniqueArray m_face_nb_nodes_in_face_mesh;
using NodeUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with node uids
using FaceUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with face uids
using EdgeUidToIndexMap = Int32UniqueArray; // use a map when no longer possible to index with edge uids
Expand All @@ -196,6 +204,7 @@ class VtkPolyhedralMeshIOService
vtkCellData* m_cell_data = nullptr;
vtkPointData* m_point_data = nullptr;
vtkCellData* m_face_data = nullptr;
vtkCellArray* m_poly_data = nullptr; // to store faces from face mesh file

void _printMeshInfos() const;

Expand All @@ -207,6 +216,7 @@ class VtkPolyhedralMeshIOService
void _checkVtkGrid() const;
void _readPlainTextVtkFaceGrid(const String& faces_filename);
void _readXmlVtkFaceGrid(const String& faces_filename);
void _readfaceNodesInFaceMesh();
};

public:
Expand Down Expand Up @@ -676,8 +686,8 @@ _createVariable(vtkDataArray* item_values, const String& variable_name, IMesh* m
void VtkPolyhedralMeshIOService::
_computeFaceVtkArcaneLidConversion(Int32Span face_vtk_to_arcane_lids, Int32Span arcane_to_vtk_lids, VtkPolyhedralMeshIOService::VtkReader& reader, IPrimaryMesh* mesh) const
{
auto face_nodes_unique_ids = reader.faceNodes();
auto face_nb_nodes = reader.faceNbNodes();
auto face_nodes_unique_ids = reader.faceNodesInFaceMesh();
auto face_nb_nodes = reader.faceNbNodesInFaceMesh();
auto current_face_index = 0;
auto current_face_index_in_face_nodes = 0;
Int32UniqueArray face_nodes_local_ids(face_nodes_unique_ids.size());
Expand All @@ -686,7 +696,7 @@ _computeFaceVtkArcaneLidConversion(Int32Span face_vtk_to_arcane_lids, Int32Span
auto current_face_nodes = face_nodes_local_ids.subConstView(current_face_index_in_face_nodes, current_face_nb_node);
current_face_index_in_face_nodes += current_face_nb_node;
Node face_first_node{ mesh->nodeFamily()->view()[current_face_nodes[0]] };
face_vtk_to_arcane_lids[current_face_index] = mesh_utils::getFaceFromNodesLocal(face_first_node, current_face_nodes).localId();
face_vtk_to_arcane_lids[current_face_index] = MeshUtils::getFaceFromNodesLocalId(face_first_node, current_face_nodes).localId();
++current_face_index;
}
auto vtk_lid = 0;
Expand Down Expand Up @@ -831,7 +841,7 @@ VtkReader(const String& filename, VtkPolyhedralTools::PrintInfoLevel print_info_
_readPlainTextVtkFaceGrid(faces_filename);
}
else{
faces_filename = m_filename + "faces.vtu";
faces_filename = m_filename + "faces.vtp";
ifile = std::ifstream{ faces_filename.localstr() };
if (ifile) {
_readXmlVtkFaceGrid(faces_filename);
Expand All @@ -842,7 +852,7 @@ VtkReader(const String& filename, VtkPolyhedralTools::PrintInfoLevel print_info_
faces_filename.split(faces_filename_and_extension,'.');

if (!ifile){
m_read_status.info_message = String::format("Information no face mesh given {0}{1} (.vtk or .vtu) to define face variables or groups on faces.",
m_read_status.info_message = String::format("Information no face mesh given {0}{1} (.vtk or .vtp) to define face variables or groups on faces.",
faces_filename_and_extension[0],
faces_filename_and_extension[1]);
}
Expand All @@ -854,6 +864,7 @@ VtkReader(const String& filename, VtkPolyhedralTools::PrintInfoLevel print_info_
}
else {
m_face_data = m_vtk_face_grid->GetCellData();
m_poly_data = m_vtk_face_grid->GetPolys();
}
}
else{
Expand Down Expand Up @@ -1306,6 +1317,28 @@ faceNbNodes()
/*---------------------------------------------------------------------------*/
/*---------------------------------------------------------------------------*/

Int64ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
faceNodesInFaceMesh()
{
if (m_face_node_uids_in_face_mesh.empty())
_readfaceNodesInFaceMesh();
return m_face_node_uids_in_face_mesh;
}

/*---------------------------------------------------------------------------*/
/*---------------------------------------------------------------------------*/

Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
faceNbNodesInFaceMesh()
{
if (m_face_nb_nodes_in_face_mesh.empty())
_readfaceNodesInFaceMesh();
return m_face_nb_nodes_in_face_mesh;
}

/*---------------------------------------------------------------------------*/
/*---------------------------------------------------------------------------*/

Int32ConstArrayView VtkPolyhedralMeshIOService::VtkReader::
edgeNbNodes()
{
Expand Down Expand Up @@ -1694,6 +1727,35 @@ _readXmlVtkFaceGrid(const String& faces_filename)
/*---------------------------------------------------------------------------*/
/*---------------------------------------------------------------------------*/

void VtkPolyhedralMeshIOService::VtkReader::
_readfaceNodesInFaceMesh()
{
m_face_nb_nodes_in_face_mesh.resize(m_poly_data->GetNumberOfCells());
m_face_node_uids_in_face_mesh.reserve(m_poly_data->GetNumberOfCells() * m_poly_data->GetMaxCellSize());
m_poly_data->InitTraversal();
vtkIdType face_nb_nodes;
vtkIdType* face_nodes;

auto face_nb_node_index = 0;
Int64UniqueArray current_face_node_uids;
Int64UniqueArray reordered_current_face_node_uids;
current_face_node_uids.reserve(m_poly_data->GetMaxCellSize());
reordered_current_face_node_uids.reserve(m_poly_data->GetMaxCellSize());
Int64UniqueArray reordered_face_node_uids(m_poly_data->GetMaxCellSize());
while (m_poly_data->GetNextCell(face_nb_nodes, face_nodes)) {
m_face_nb_nodes_in_face_mesh[face_nb_node_index] = face_nb_nodes;
ConstArrayView<vtkIdType> face_nodes_view(face_nb_nodes, face_nodes);
current_face_node_uids.resize(face_nb_nodes);
reordered_current_face_node_uids.resize(face_nb_nodes);
std::copy(face_nodes_view.begin(), face_nodes_view.end(), current_face_node_uids.begin());
MeshUtils::reorderNodesOfFace(current_face_node_uids, reordered_current_face_node_uids);
std::copy(reordered_current_face_node_uids.begin(), reordered_current_face_node_uids.end(), std::back_inserter(m_face_node_uids_in_face_mesh));
++face_nb_node_index;
}
}
/*---------------------------------------------------------------------------*/
/*---------------------------------------------------------------------------*/

void VtkPolyhedralMeshIOService::VtkReader::
_checkVtkGrid() const
{
Expand Down
2 changes: 1 addition & 1 deletion arcane/src/arcane/tests/CMakeLists.txt
Original file line number Diff line number Diff line change
Expand Up @@ -263,7 +263,7 @@ arcane_copy_mesh_direct(cross_a_2x1x1.vtkfaces.vtk)
arcane_copy_mesh_direct(fan_2x1x3.vtk)
arcane_copy_mesh_direct(fan_2x1x3.vtkfaces.vtk)
arcane_copy_mesh_direct(cross_a_2x1x1.vtu)
arcane_copy_mesh_direct(cross_a_2x1x1.vtufaces.vtu)
arcane_copy_mesh_direct(cross_a_2x1x1.vtufaces.vtp)
arcane_copy_mesh_direct(plancher.msh)
arcane_copy_mesh_direct(plancher.binary.msh)
arcane_copy_mesh_direct(hex_tetra_pyramics.msh)
Expand Down
Loading